Geant4 10.7.0
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 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)
 
virtual ~G4VITProcess ()
 
 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 ()
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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 void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- 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 49 of file G4DNAMolecularDissociation.cc.

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

◆ G4DNAMolecularDissociation() [2/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( )
delete

◆ G4DNAMolecularDissociation() [3/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4DNAMolecularDissociation right)
delete

◆ ~G4DNAMolecularDissociation()

G4DNAMolecularDissociation::~G4DNAMolecularDissociation ( )
overridedefault

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 314 of file G4DNAMolecularDissociation.cc.

316{
319 return DecayIt(track, step);
320}
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

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 324 of file G4DNAMolecularDissociation.cc.

326{
327 if (fDecayAtFixedTime)
328 {
329 return GetMeanLifeTime(track, condition);
330 }
331
333}
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 115 of file G4DNAMolecularDissociation.cc.

117{
119 auto pMotherMolecule = GetMolecule(track);
120 auto pMotherMoleculeDefinition = pMotherMolecule->GetDefinition();
121
122 if (pMotherMoleculeDefinition->GetDecayTable())
123 {
124 const auto pDissociationChannels = pMotherMolecule->GetDissociationChannels();
125
126 if (pDissociationChannels == nullptr)
127 {
128 G4ExceptionDescription exceptionDescription;
129 pMotherMolecule->PrintState();
130 exceptionDescription << "No decay channel was found for the molecule : "
131 << pMotherMolecule->GetName() << G4endl;
132 G4Exception("G4DNAMolecularDissociation::DecayIt",
133 "G4DNAMolecularDissociation::NoDecayChannel",
135 exceptionDescription);
136 return &aParticleChange;
137 }
138
139 auto decayVectorSize = pDissociationChannels->size();
140 G4double RdmValue = G4UniformRand();
141
142 const G4MolecularDissociationChannel* pDecayChannel = nullptr;
143 size_t i = 0;
144 do
145 {
146 pDecayChannel = (*pDissociationChannels)[i];
147 if (RdmValue < pDecayChannel->GetProbability())
148 {
149 break;
150 }
151 RdmValue -= pDecayChannel->GetProbability();
152 i++;
153 } while (i < decayVectorSize);
154
155 G4double decayEnergy = pDecayChannel->GetEnergy();
156 auto nbProducts = pDecayChannel->GetNbProducts();
157
158 if (decayEnergy > 0.)
159 {
161 }
162
163 if (nbProducts)
164 {
165 std::vector<G4ThreeVector> productsDisplacement(nbProducts);
166 G4ThreeVector motherMoleculeDisplacement;
167
168 auto it = fDisplacementMap.find(pMotherMoleculeDefinition);
169
170 if (it != fDisplacementMap.end())
171 {
172 auto pDisplacer = it->second.get();
173 productsDisplacement = pDisplacer->GetProductsDisplacement(pDecayChannel);
174 motherMoleculeDisplacement =
175 pDisplacer->GetMotherMoleculeDisplacement(pDecayChannel);
176 }
177 else
178 {
180 errMsg << "No G4MolecularDecayProcess::theDecayDisplacementMap["
181 << pMotherMolecule->GetName() + "]";
182 G4Exception("G4MolecularDecayProcess::DecayIt",
183 "DNAMolecularDecay001",
185 errMsg);
186 }
187
189
190#ifdef G4VERBOSE
191 if (fVerbose)
192 {
193 G4cout << "Decay Process : " << pMotherMolecule->GetName()
194 << " (trackID :" << track.GetTrackID() << ") "
195 << pDecayChannel->GetName() << G4endl;
196 }
197#endif
198
200
201 for (G4int j = 0; j < nbProducts; j++)
202 {
203 auto pProduct = new G4Molecule(pDecayChannel->GetProduct(j));
204
205 G4ThreeVector displacement = motherMoleculeDisplacement + productsDisplacement[j];
206 double mag_displacement = displacement.mag();
207 G4ThreeVector displacement_direction = displacement / (mag_displacement + 1e-30);
208
209 double prNewSafety = DBL_MAX;
210
211 //double step =
212 pNavigator->CheckNextStep(track.GetPosition(),
213 displacement_direction,
214 mag_displacement,
215 prNewSafety);
216
217 //if(prNewSafety < mag_displacement || step < mag_displacement)
218 mag_displacement = std::min(prNewSafety * 0.8, mag_displacement);
219
220 G4ThreeVector product_pos = track.GetPosition()
221 + displacement_direction * mag_displacement;
222
223 const G4AffineTransform& transform = pNavigator->GetGlobalToLocalTransform();
224
225 G4ThreeVector localPoint = transform.TransformPoint(product_pos); //track.GetPosition());
226
227 if (track.GetTouchable()->GetSolid()->Inside(localPoint) !=
228 EInside::kInside)
229 {
231 ED << "The decayed product is outside of the volume : "
232 << track.GetTouchable()->GetVolume()->GetName() << G4endl;
233 G4Exception("G4DNAMolecularDissociation::DecayIt()",
234 "OUTSIDE_OF_MOTHER_VOLUME",
235 JustWarning, ED);
236 }
237
238 auto pSecondary = pProduct->BuildTrack(track.GetGlobalTime(), product_pos);
239
240 pSecondary->SetTrackStatus(fAlive);
241#ifdef G4VERBOSE
242 if (fVerbose)
243 {
244 G4cout << "Product : " << pProduct->GetName() << G4endl;
245 }
246#endif
247 // add the secondary track in the List
248 aParticleChange.G4VParticleChange::AddSecondary(pSecondary);
249 }
250#ifdef G4VERBOSE
251 if (fVerbose)
252 {
253 G4cout << "-------------" << G4endl;
254 }
255#endif
256 }
257 //DEBUG
258 else if (fVerbose && decayEnergy)
259 {
260 G4cout << "No products for this channel" << G4endl;
261 G4cout << "-------------" << G4endl;
262 }
263 /*
264 else if(!decayEnergy && !nbProducts)
265 {
266 G4ExceptionDescription errMsg;
267 errMsg << "There is no products and no energy specified in the molecular "
268 "dissociation channel";
269 G4Exception("G4MolecularDissociationProcess::DecayIt",
270 "DNAMolecularDissociation002",
271 FatalErrorInArgument,
272 errMsg);
273 }
274 */
275 }
276
278
279 return &aParticleChange;
280}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
@ 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
virtual void Initialize(const G4Track &)
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 G4VSolid * GetSolid(G4int depth=0) const
Definition: G4VTouchable.cc:48
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41
#define DBL_MAX
Definition: templates.hh:62

Referenced by AtRestDoIt().

◆ GetDisplacer()

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

Definition at line 291 of file G4DNAMolecularDissociation.cc.

292{
293 return fDisplacementMap[pSpecies].get();
294}

◆ GetMeanFreePath()

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

Implements G4VITRestDiscreteProcess.

Definition at line 345 of file G4DNAMolecularDissociation.cc.

348{
349 return 0;
350}

◆ GetMeanLifeTime()

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

Implements G4VITRestDiscreteProcess.

Definition at line 106 of file G4DNAMolecularDissociation.cc.

108{
109 G4double output = GetMolecule(track)->GetDecayTime() - track.GetProperTime();
110 return output > 0. ? output : 0.;
111}
G4double GetDecayTime() const
Definition: G4Molecule.cc:474
G4double GetProperTime() const

Referenced by AtRestGetPhysicalInteractionLength().

◆ IsApplicable()

G4bool G4DNAMolecularDissociation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 82 of file G4DNAMolecularDissociation.cc.

84{
85 if (aParticleType.GetParticleType() == "Molecule")
86 {
87#ifdef G4VERBOSE
88
89 if (fVerbose > 1)
90 {
91 G4cout << "G4MolecularDissociation::IsApplicable(";
92 G4cout << aParticleType.GetParticleName() << ",";
93 G4cout << aParticleType.GetParticleType() << ")" << G4endl;
94 }
95#endif
96 return (true);
97 }
98 else
99 {
100 return false;
101 }
102}
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

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 337 of file G4DNAMolecularDissociation.cc.

339{
340 return AtRestDoIt(track, step);
341}
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &step) override

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 298 of file G4DNAMolecularDissociation.cc.

301{
302 return 0; //-1*(-log(1-G4UniformRand())*100*1e-15*s);
303}

◆ SetDisplacer()

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

Definition at line 284 of file G4DNAMolecularDissociation.cc.

285{
286 fDisplacementMap.emplace(pSpecies, std::unique_ptr<Displacer>(pDisplacer));
287}

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

◆ SetVerbose()

void G4DNAMolecularDissociation::SetVerbose ( G4int  verbose)

Definition at line 307 of file G4DNAMolecularDissociation.cc.

308{
309 fVerbose = verbose;
310}

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