Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4DNAMolecularDissociation Class Reference

#include <G4DNAMolecularDissociation.hh>

Inheritance diagram for G4DNAMolecularDissociation:
G4VITRestProcess G4VITProcess G4VProcess

Public Member Functions

 G4DNAMolecularDissociation (const G4String &processName="DNAMolecularDecay", G4ProcessType type=fDecay)
 
virtual ~G4DNAMolecularDissociation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &step)
 
void SetVerbose (G4int)
 
void SetDecayDisplacer (const G4ParticleDefinition *, G4VMolecularDecayDisplacer *)
 
G4VMolecularDecayDisplacerGetDecayDisplacer (const G4ParticleDefinition *)
 
- Public Member Functions inherited from G4VITRestProcess
 G4VITRestProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITRestProcess (const G4VITRestProcess &)
 
virtual ~G4VITRestProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
virtual ~G4VITProcess ()
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4int operator== (const G4VITProcess &right) const
 
G4int operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4ProcessState_Lock *aProcInfo)
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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 &)
 
virtual G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
- Protected Member Functions inherited from G4VITRestProcess
 G4VITRestProcess ()
 
G4VITRestProcessoperator= (const G4VITRestProcess &right)
 
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
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
G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

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 56 of file G4DNAMolecularDissociation.hh.

Constructor & Destructor Documentation

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

Definition at line 48 of file G4DNAMolecularDissociation.cc.

References G4VProcess::aParticleChange, G4VProcess::enableAlongStepDoIt, G4VProcess::enableAtRestDoIt, G4VProcess::enablePostStepDoIt, G4cout, G4DNAMolecularDissociation(), G4endl, G4VProcess::pParticleChange, G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.

Referenced by G4DNAMolecularDissociation().

49  : G4VITRestProcess(processName, type)
50 {
51  // set Process Sub Type
52  SetProcessSubType(59); // DNA sub-type
53  enableAlongStepDoIt = false;
54  enablePostStepDoIt = false;
55  enableAtRestDoIt=true;
56 
57  fVerbose = 0 ;
58 
59 #ifdef G4VERBOSE
60  if (verboseLevel>1)
61  {
62  G4cout << "G4MolecularDecayProcess constructor " << " Name:" << processName << G4endl;
63  }
64 #endif
65 
67 
68  fDecayAtFixedTime = true ;
69 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4DNAMolecularDissociation::~G4DNAMolecularDissociation ( )
virtual

Definition at line 71 of file G4DNAMolecularDissociation.cc.

72 {
73  DecayDisplacementMap::iterator it = fDecayDisplacementMap.begin();
74 
75  for( ; it != fDecayDisplacementMap.end() ; it++)
76  {
77  if(it->second)
78  {
79  delete it->second;
80  it->second = 0;
81  }
82  }
83  fDecayDisplacementMap.clear();
84 }

Member Function Documentation

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

Reimplemented from G4VITRestProcess.

Definition at line 120 of file G4DNAMolecularDissociation.hh.

References G4VITProcess::ClearInteractionTimeLeft(), G4VITProcess::ClearNumberOfInteractionLengthLeft(), and DecayIt().

124 {
127  return DecayIt(track, step);
128 }
virtual void ClearNumberOfInteractionLengthLeft()
virtual void ClearInteractionTimeLeft()
virtual G4VParticleChange * DecayIt(const G4Track &, const G4Step &)
G4double G4DNAMolecularDissociation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
inlinevirtual

Reimplemented from G4VITRestProcess.

Definition at line 108 of file G4DNAMolecularDissociation.hh.

References G4VITRestProcess::AtRestGetPhysicalInteractionLength(), and GetMeanLifeTime().

111 {
112  if(fDecayAtFixedTime)
113  {
114  return GetMeanLifeTime(track, condition);
115  }
116 
118 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
G4VParticleChange * G4DNAMolecularDissociation::DecayIt ( const G4Track track,
const G4Step  
)
protectedvirtual

Definition at line 121 of file G4DNAMolecularDissociation.cc.

References G4VProcess::aParticleChange, G4Molecule::BuildTrack(), G4ElectronOccupancy::DumpInfo(), fAlive, FatalErrorInArgument, FatalException, fStopAndKill, G4cout, G4endl, G4Exception(), G4UniformRand, G4Molecule::GetDefinition(), G4Molecule::GetElectronOccupancy(), G4MolecularDecayChannel::GetEnergy(), GetMolecule(), G4MolecularDecayChannel::GetName(), G4Molecule::GetName(), G4MolecularDecayChannel::GetNbProducts(), G4Track::GetPosition(), G4MolecularDecayChannel::GetProbability(), G4MolecularDecayChannel::GetProduct(), G4VMolecularDecayDisplacer::GetProductsDisplacement(), G4Track::GetTrackID(), G4ParticleChange::Initialize(), G4ITManager< T >::Instance(), python.hepunit::picosecond, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), and G4VParticleChange::SetNumberOfSecondaries().

Referenced by AtRestDoIt().

124 {
125  // DEBUG
126  // G4cout << "Is calling G4MolecularDecayProcess::DecayIt" << G4endl;
127 
129  const G4Molecule * theMotherMolecule = GetMolecule(track);
130  const G4MoleculeDefinition* moleculeDefinition = theMotherMolecule->GetDefinition();
131 
132  // DEBUG
133  // G4cout <<"Calling G4MolecularDecayProcess::DecayIt"<<G4endl;
134  // G4cout << "The mother molecule state : " << G4endl;
135  // theMotherMolecule -> PrintState();
136 
137  if(moleculeDefinition-> GetDecayTable())
138  {
139  const vector<const G4MolecularDecayChannel*>* DecayVector =
140  (theMotherMolecule -> GetDecayChannel());
141 
142  if(DecayVector == 0)
143  {
144  G4ExceptionDescription exceptionDescription;
145  theMotherMolecule->GetElectronOccupancy()->DumpInfo();
146  exceptionDescription << "No decay channel was found for the molecule : " << theMotherMolecule-> GetName() << G4endl;
147  G4Exception("G4DNAMolecularDissociation::DecayIt", "G4DNAMolecularDissociation::NoDecayChannel",FatalException,exceptionDescription);
148  return &aParticleChange;
149  }
150 
151  G4int DecayVectorSize = DecayVector-> size();
152  // DEBUG
153  // G4cout<< "Number of decay channels : " << DecayVectorSize<<G4endl;
154  G4double RdmValue = G4UniformRand();
155 
156  const G4MolecularDecayChannel* decayChannel(0);
157  G4int i=0;
158  do
159  {
160  decayChannel = (*DecayVector)[i];
161  if(RdmValue < decayChannel->GetProbability()) break;
162  RdmValue -= decayChannel->GetProbability();
163  i++;
164  }
165  while(i< DecayVectorSize);
166 
167  // DEBUG
168  // G4cout<< "Selected Decay channel : " << decayChannel->GetName()<<G4endl;
169 
170  G4double decayEnergy = decayChannel->GetEnergy();
171  G4int nbProducts = decayChannel->GetNbProducts();
172 
173  if(decayEnergy)
174  {
175  // DEBUG
176  // G4cout<<"Deposit energy :" <<decayChannel->GetEnergy()/eV << " eV"<<G4endl;
177 
178  aParticleChange.ProposeLocalEnergyDeposit(decayChannel->GetEnergy());
179  }
180 
181  if(nbProducts)
182  {
183 
184  // DEBUG
185  // G4cout<<"Number of products :" <<nbProducts<<G4endl;
186 
187  vector<G4ThreeVector> ProductsDisplacement(nbProducts);
188  G4ThreeVector theMotherMoleculeDisplacement;
189 
190  DecayDisplacementMap::iterator it = fDecayDisplacementMap.find(moleculeDefinition);
191 
192  if(it!=fDecayDisplacementMap.end())
193  {
194  G4VMolecularDecayDisplacer* displacer = it->second;
195  ProductsDisplacement = displacer->GetProductsDisplacement(decayChannel);
196  theMotherMoleculeDisplacement = displacer-> GetMotherMoleculeDisplacement(decayChannel);
197  }
198  else
199  {
200  G4ExceptionDescription errMsg;
201  errMsg << "No G4MolecularDecayProcess::theDecayDisplacementMap["
202  << theMotherMolecule->GetName() +"]" ;
203  G4Exception("G4MolecularDecayProcess::DecayIt","DNAMolecularDecay001",FatalErrorInArgument, errMsg);
204  }
205 
207 
208 #ifdef G4VERBOSE
209  if(fVerbose)
210  {
211  G4cout<<"Decay Process : "
212  << theMotherMolecule->GetName()
213  << " (trackID :" << track.GetTrackID() << ") "
214  << decayChannel->GetName()
215  << G4endl;
216  }
217 #endif
218 
219  for (G4int j=0; j<nbProducts ; j++)
220  {
221  G4Molecule* product = new G4Molecule(*decayChannel->GetProduct(j));
222 
223  // create a new track object
224  // Be carefull as this processes is dedicated to be used in the DNA module
225  // The molecular decay will happen one picosecond after the start of the simulation.
226  // This may be seen as a bug and will be hopefully improve in the next releases
227  G4Track* secondary = product->BuildTrack(picosecond,track.GetPosition()
228  + theMotherMoleculeDisplacement + ProductsDisplacement[j]);
229 
230  secondary-> SetTrackStatus(fAlive);
231 #ifdef G4VERBOSE
232  if(fVerbose)
233  {
234  G4cout<<"Product : "<< product->GetName()<<G4endl;
235  }
236 #endif
237  // add the secondary track in the List
238  aParticleChange.G4VParticleChange::AddSecondary(secondary);
239  G4ITManager<G4Molecule>::Instance()->Push(secondary);
240  }
241 #ifdef G4VERBOSE
242  if(fVerbose)
243  G4cout<<"-------------"<<G4endl;
244 #endif
245  }
246  // DEBUG
247  // else if(decayEnergy)
248  // {
249  // G4cout << "No products for this channel" << G4endl ;
250  // G4cout<<"-------------"<<G4endl;
251  // }
252  else if(!decayEnergy && !nbProducts)
253  {
254  G4ExceptionDescription errMsg;
255  errMsg << "There is no products and no energy specified in the molecular decay channel";
256  G4Exception("G4MolecularDecayProcess::DecayIt","DNAMolecularDecay002",FatalErrorInArgument, errMsg);
257  }
258  }
259 
261 
262  return &aParticleChange;
263 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4ThreeVector & GetPosition() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4String & GetName() const
Definition: G4Molecule.cc:243
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
static G4ITManager< T > * Instance()
virtual std::vector< G4ThreeVector > GetProductsDisplacement(const G4MolecularDecayChannel *) const =0
G4int GetTrackID() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
const G4MoleculeDefinition * GetDefinition() const
Definition: G4Molecule.cc:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4Track * BuildTrack(G4double globalTime, const G4ThreeVector &Position)
Definition: G4Molecule.cc:263
virtual void Initialize(const G4Track &)
int picosecond
Definition: hepunit.py:86
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4ElectronOccupancy * GetElectronOccupancy() const
Definition: G4Molecule.cc:374
G4VMolecularDecayDisplacer * G4DNAMolecularDissociation::GetDecayDisplacer ( const G4ParticleDefinition molDef)

Definition at line 270 of file G4DNAMolecularDissociation.cc.

271 {
272  return fDecayDisplacementMap[molDef] ;
273 }
G4double G4DNAMolecularDissociation::GetMeanLifeTime ( const G4Track track,
G4ForceCondition  
)
protectedvirtual

Implements G4VITRestProcess.

Definition at line 114 of file G4DNAMolecularDissociation.cc.

References GetMolecule(), and G4Track::GetProperTime().

Referenced by AtRestGetPhysicalInteractionLength().

116 {
117  G4double output = GetMolecule(track)-> GetDecayTime() - track.GetProperTime() ;
118  return (output > 0 ? output : 0 );
119 }
G4double GetProperTime() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
double G4double
Definition: G4Types.hh:76
G4bool G4DNAMolecularDissociation::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 94 of file G4DNAMolecularDissociation.cc.

References G4cout, G4endl, G4ParticleDefinition::GetParticleName(), and G4ParticleDefinition::GetParticleType().

95 {
96  if(aParticleType.GetParticleType()=="Molecule")
97  {
98 #ifdef G4VERBOSE
99  if(fVerbose>1)
100  {
101  G4cout<<"G4MolecularDecay::IsApplicable(";
102  G4cout<<aParticleType.GetParticleName()<<",";
103  G4cout<<aParticleType.GetParticleType()<<")"<<G4endl;
104  }
105 #endif
106  return(true);
107  }
108  else
109  {
110  return false;
111  }
112 }
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleType() const
#define G4endl
Definition: G4ios.hh:61
void G4DNAMolecularDissociation::SetDecayDisplacer ( const G4ParticleDefinition molDef,
G4VMolecularDecayDisplacer aDisplacer 
)

Definition at line 265 of file G4DNAMolecularDissociation.cc.

266 {
267  fDecayDisplacementMap[molDef] = aDisplacer;
268 }
void G4DNAMolecularDissociation::SetVerbose ( G4int  verbose)
inline

Definition at line 103 of file G4DNAMolecularDissociation.hh.

104 {
105  fVerbose = verbose ;
106 }

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