Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMolecularDissociation Class Reference

#include <G4DNAMolecularDissociation.hh>

+ Inheritance diagram for G4DNAMolecularDissociation:

Public Types

using Species = const G4MoleculeDefinition
 
using Displacer = G4VMolecularDissociationDisplacer
 

Public Member Functions

 G4DNAMolecularDissociation (const G4String &processName="DNAMolecularDecay", G4ProcessType type=fDecay)
 
 G4DNAMolecularDissociation ()=delete
 
 G4DNAMolecularDissociation (const G4DNAMolecularDissociation &right)=delete
 
G4DNAMolecularDissociationoperator= (const G4DNAMolecularDissociation &right)=delete
 
 ~G4DNAMolecularDissociation () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &step) override
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step) override
 
void SetVerbose (G4int)
 
void SetUserBrownianAction (G4VUserBrownianAction *pBrownianAction)
 
void SetDisplacer (Species *, Displacer *)
 
DisplacerGetDisplacer (Species *)
 
- Public Member Functions inherited from G4VITRestDiscreteProcess
 G4VITRestDiscreteProcess ()=delete
 
 G4VITRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITRestDiscreteProcess (const G4VITRestDiscreteProcess &)=delete
 
G4VITRestDiscreteProcessoperator= (const G4VITRestDiscreteProcess &right)=delete
 
 ~G4VITRestDiscreteProcess () override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
 ~G4VITProcess () override
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
void StartTracking (G4Track *) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double GetInteractionTimeLeft ()
 
void ResetNumberOfInteractionLengthLeft () override
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &, const G4Step &)
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *) override
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *) override
 
- Protected Member Functions inherited from G4VITRestDiscreteProcess
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
template<typename T >
T * GetState ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VITProcess
G4shared_ptr< G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

G4DNAMolecularDissociation should be called only for molecules. It will dissociate the molecules using the decay associated to this molecule and if a displacement scheme has been registered, it will place the products to the expected position.

Definition at line 63 of file G4DNAMolecularDissociation.hh.

Member Typedef Documentation

◆ Displacer

◆ Species

Constructor & Destructor Documentation

◆ G4DNAMolecularDissociation() [1/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4String & processName = "DNAMolecularDecay",
G4ProcessType type = fDecay )

Definition at line 50 of file G4DNAMolecularDissociation.cc.

53 : G4VITRestDiscreteProcess(processName, type)
54{
55 // set Process Sub Type
57 enableAlongStepDoIt = false;
58 enablePostStepDoIt = true;
59 enableAtRestDoIt = true;
60
61 fVerbose = 0;
62
63#ifdef G4VERBOSE
64 if (verboseLevel > 1)
65 {
66 G4cout << "G4MolecularDissociationProcess constructor " << " Name:"
67 << processName << G4endl;
68 }
69#endif
70
72
73 fDecayAtFixedTime = true;
74 fProposesTimeStep = true;
75}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4bool fProposesTimeStep
G4ParticleChange aParticleChange
G4bool enableAtRestDoIt
G4int verboseLevel
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange

◆ G4DNAMolecularDissociation() [2/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( )
delete

◆ G4DNAMolecularDissociation() [3/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4DNAMolecularDissociation & right)
delete

◆ ~G4DNAMolecularDissociation()

G4DNAMolecularDissociation::~G4DNAMolecularDissociation ( )
override

Definition at line 79 of file G4DNAMolecularDissociation.cc.

80{
81 delete fpBrownianAction;
82}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4DNAMolecularDissociation::AtRestDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 334 of file G4DNAMolecularDissociation.cc.

336{
339 return DecayIt(track, step);
340}
virtual G4VParticleChange * DecayIt(const G4Track &, const G4Step &)
virtual void ClearInteractionTimeLeft()
virtual void ClearNumberOfInteractionLengthLeft()

Referenced by PostStepDoIt().

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 344 of file G4DNAMolecularDissociation.cc.

346{
347 if (fDecayAtFixedTime)
348 {
349 return GetMeanLifeTime(track, condition);
350 }
351
353}
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *) override
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override

◆ DecayIt()

G4VParticleChange * G4DNAMolecularDissociation::DecayIt ( const G4Track & track,
const G4Step &  )
protectedvirtual

Definition at line 117 of file G4DNAMolecularDissociation.cc.

119{
121 auto pMotherMolecule = GetMolecule(track);
122 auto pMotherMoleculeDefinition = pMotherMolecule->GetDefinition();
123
124 if (pMotherMoleculeDefinition->GetDecayTable() != nullptr)
125 {
126 const auto pDissociationChannels = pMotherMolecule->GetDissociationChannels();
127
128 if (pDissociationChannels == nullptr)
129 {
130 G4ExceptionDescription exceptionDescription;
131 pMotherMolecule->PrintState();
132 exceptionDescription << "No decay channel was found for the molecule : "
133 << pMotherMolecule->GetName() << G4endl;
134 G4Exception("G4DNAMolecularDissociation::DecayIt",
135 "G4DNAMolecularDissociation::NoDecayChannel",
137 exceptionDescription);
138 return &aParticleChange;
139 }
140
141 auto decayVectorSize = pDissociationChannels->size();
142 G4double RdmValue = G4UniformRand();
143
144 const G4MolecularDissociationChannel* pDecayChannel = nullptr;
145 size_t i = 0;
146 do
147 {
148 pDecayChannel = (*pDissociationChannels)[i];
149 if (RdmValue < pDecayChannel->GetProbability())
150 {
151 break;
152 }
153 RdmValue -= pDecayChannel->GetProbability();
154 i++;
155 } while (i < decayVectorSize);
156
157 G4double decayEnergy = pDecayChannel->GetEnergy();
158 auto nbProducts = pDecayChannel->GetNbProducts();
159
160 if (decayEnergy > 0.)
161 {
163 }
164
165 if (nbProducts != 0)
166 {
167 std::vector<G4ThreeVector> productsDisplacement(nbProducts);
168 G4ThreeVector motherMoleculeDisplacement;
169
170 auto it = fDisplacementMap.find(pMotherMoleculeDefinition);
171
172 if (it != fDisplacementMap.end())
173 {
174 auto pDisplacer = it->second.get();
175 productsDisplacement = pDisplacer->GetProductsDisplacement(pDecayChannel);
176 motherMoleculeDisplacement =
177 pDisplacer->GetMotherMoleculeDisplacement(pDecayChannel);
178 }
179 else
180 {
182 errMsg << "No G4MolecularDecayProcess::theDecayDisplacementMap["
183 << pMotherMolecule->GetName() + "]";
184 G4Exception("G4MolecularDecayProcess::DecayIt",
185 "DNAMolecularDecay001",
187 errMsg);
188 }
189
191
192#ifdef G4VERBOSE
193 if (fVerbose != 0)
194 {
195 G4cout << "Decay Process : " << pMotherMolecule->GetName()
196 << " (trackID :" << track.GetTrackID() << ") "
197 << pDecayChannel->GetName() << G4endl;
198 }
199#endif
200
202
203 for (G4int j = 0; j < nbProducts; j++)
204 {
205 auto pProduct = new G4Molecule(pDecayChannel->GetProduct(j));
206
207 G4ThreeVector displacement = motherMoleculeDisplacement + productsDisplacement[j];
208 double mag_displacement = displacement.mag();
209 G4ThreeVector displacement_direction = displacement / (mag_displacement + 1e-30);
210
211 double prNewSafety = DBL_MAX;
212
213 //double step =
214 pNavigator->CheckNextStep(track.GetPosition(),
215 displacement_direction,
216 mag_displacement,
217 prNewSafety);
218
219 //if(prNewSafety < mag_displacement || step < mag_displacement)
220 mag_displacement = std::min(prNewSafety * 0.8, mag_displacement);
221
222 G4ThreeVector product_pos = track.GetPosition()
223 + displacement_direction * mag_displacement;
224
225 //Hoang: force changing position track::
226 if(fpBrownianAction != nullptr)
227 {
228 fpBrownianAction->Transport(product_pos);
229 }
230 //Hoang: force changing position track
231
232 const G4AffineTransform& transform = pNavigator->GetGlobalToLocalTransform();
233
234 G4ThreeVector localPoint = transform.TransformPoint(product_pos); //track.GetPosition());
235 // warning if the decayed product is outside of the volume and
236 // the mother volume has no water material (the decayed product
237 // is outside of the world volume will be killed in the next step)
238 if (track.GetTouchable()->GetSolid()->Inside(localPoint) !=
240 {
241 auto WaterMaterial = G4Material::GetMaterial("G4_WATER");
242 auto Motherlogic = track.GetTouchable()->GetVolume()->
243 GetMotherLogical();
244 if (Motherlogic != nullptr
245 && Motherlogic->GetMaterial() != WaterMaterial)
246 {
248 ED << "The decayed product is outside of the volume : "
249 << track.GetTouchable()->GetVolume()->GetName()
250 << " with material : "<< Motherlogic->GetMaterial()
251 ->GetName()<< G4endl;
252 G4Exception("G4DNAMolecularDissociation::DecayIt()",
253 "OUTSIDE_OF_MOTHER_VOLUME",
254 JustWarning, ED);
255 }
256 }
257
258 auto pSecondary = pProduct->BuildTrack(track.GetGlobalTime(), product_pos);
259
260 pSecondary->SetTrackStatus(fAlive);
261#ifdef G4VERBOSE
262 if (fVerbose != 0)
263 {
264 G4cout << "Product : " << pProduct->GetName() << G4endl;
265 }
266#endif
267 // add the secondary track in the List
268 aParticleChange.G4VParticleChange::AddSecondary(pSecondary);
269 }
270#ifdef G4VERBOSE
271 if (fVerbose != 0)
272 {
273 G4cout << "-------------" << G4endl;
274 }
275#endif
276 }
277 //DEBUG
278 else if ((fVerbose != 0) && (decayEnergy != 0.0))
279 {
280 G4cout << "No products for this channel" << G4endl;
281 G4cout << "-------------" << G4endl;
282 }
283 /*
284 else if(!decayEnergy && !nbProducts)
285 {
286 G4ExceptionDescription errMsg;
287 errMsg << "There is no products and no energy specified in the molecular "
288 "dissociation channel";
289 G4Exception("G4MolecularDissociationProcess::DecayIt",
290 "DNAMolecularDissociation002",
291 FatalErrorInArgument,
292 errMsg);
293 }
294 */
295 }
296
298
299 return &aParticleChange;
300}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
@ fAlive
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
double mag() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
static G4ITTransportationManager * GetTransportationManager()
G4ITNavigator * GetNavigatorForTracking() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
void Initialize(const G4Track &) override
virtual G4VSolid * GetSolid(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4VTouchable * GetTouchable() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual void Transport(G4ThreeVector &, G4Track *pTrack=nullptr)=0
@ kInside
Definition geomdefs.hh:70
#define DBL_MAX
Definition templates.hh:62

Referenced by AtRestDoIt().

◆ GetDisplacer()

G4DNAMolecularDissociation::Displacer * G4DNAMolecularDissociation::GetDisplacer ( Species * pSpecies)

Definition at line 311 of file G4DNAMolecularDissociation.cc.

312{
313 return fDisplacementMap[pSpecies].get();
314}

◆ GetMeanFreePath()

G4double G4DNAMolecularDissociation::GetMeanFreePath ( const G4Track & ,
G4double ,
G4ForceCondition *  )
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 365 of file G4DNAMolecularDissociation.cc.

368{
369 return 0;
370}

◆ GetMeanLifeTime()

G4double G4DNAMolecularDissociation::GetMeanLifeTime ( const G4Track & track,
G4ForceCondition *  )
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 108 of file G4DNAMolecularDissociation.cc.

110{
111 G4double output = GetMolecule(track)->GetDecayTime() - track.GetProperTime();
112 return output > 0. ? output : 0.;
113}
G4double GetDecayTime() const
G4double GetProperTime() const

Referenced by AtRestGetPhysicalInteractionLength().

◆ IsApplicable()

G4bool G4DNAMolecularDissociation::IsApplicable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4DNAMolecularDissociation.cc.

88{
89 if (aParticleType.GetParticleType() == "Molecule")
90 {
91#ifdef G4VERBOSE
92
93 if (fVerbose > 1)
94 {
95 G4cout << "G4MolecularDissociation::IsApplicable(";
96 G4cout << aParticleType.GetParticleName() << ",";
97 G4cout << aParticleType.GetParticleType() << ")" << G4endl;
98 }
99#endif
100 return (true);
101 }
102
103 return false;
104}
const G4String & GetParticleType() const
const G4String & GetParticleName() const

◆ operator=()

G4DNAMolecularDissociation & G4DNAMolecularDissociation::operator= ( const G4DNAMolecularDissociation & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4DNAMolecularDissociation::PostStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 357 of file G4DNAMolecularDissociation.cc.

359{
360 return AtRestDoIt(track, step);
361}
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &step) override

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 318 of file G4DNAMolecularDissociation.cc.

321{
322 return 0; //-1*(-log(1-G4UniformRand())*100*1e-15*s);
323}

◆ SetDisplacer()

void G4DNAMolecularDissociation::SetDisplacer ( Species * pSpecies,
Displacer * pDisplacer )

Definition at line 304 of file G4DNAMolecularDissociation.cc.

305{
306 fDisplacementMap.emplace(pSpecies, std::unique_ptr<Displacer>(pDisplacer));
307}

Referenced by G4EmDNAChemistry::ConstructProcess(), G4EmDNAChemistry_option1::ConstructProcess(), and G4EmDNAChemistry_option2::ConstructProcess().

◆ SetUserBrownianAction()

void G4DNAMolecularDissociation::SetUserBrownianAction ( G4VUserBrownianAction * pBrownianAction)
inline

Definition at line 91 of file G4DNAMolecularDissociation.hh.

92 {
93 fpBrownianAction = pBrownianAction;
94 }

◆ SetVerbose()

void G4DNAMolecularDissociation::SetVerbose ( G4int verbose)

Definition at line 327 of file G4DNAMolecularDissociation.cc.

328{
329 fVerbose = verbose;
330}

The documentation for this class was generated from the following files: