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

#include <G4DNAMoleculeEncounterStepper.hh>

Inheritance diagram for G4DNAMoleculeEncounterStepper:
G4VITTimeStepper

Public Member Functions

 G4DNAMoleculeEncounterStepper ()
 
virtual ~G4DNAMoleculeEncounterStepper ()
 
 G4DNAMoleculeEncounterStepper (const G4DNAMoleculeEncounterStepper &)
 
virtual void Prepare ()
 
virtual G4double CalculateStep (const G4Track &, const G4double &)
 
void SetReactionModel (G4VDNAReactionModel *)
 
G4VDNAReactionModelGetReactionModel ()
 
void SetVerbose (int)
 
- Public Member Functions inherited from G4VITTimeStepper
 G4VITTimeStepper ()
 
virtual ~G4VITTimeStepper ()
 
 G4VITTimeStepper (const G4VITTimeStepper &)
 
G4VITTimeStepperoperator= (const G4VITTimeStepper &other)
 
virtual void Initialize ()
 
G4TrackVectorHandle GetReactants ()
 
virtual void ResetReactants ()
 
G4double GetSampledMinTimeStep ()
 
void SetReactionTable (const G4ITReactionTable *)
 
const G4ITReactionTableGetReactionTable ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITTimeStepper
static void SetTimes (const G4double &, const G4double &)
 
- Protected Attributes inherited from G4VITTimeStepper
G4double fSampledMinTimeStep
 
G4TrackVectorHandle fReactants
 
const G4ITReactionTablefpReactionTable
 
- Static Protected Attributes inherited from G4VITTimeStepper
static G4ThreadLocal G4double fCurrentGlobalTime = -1
 
static G4ThreadLocal G4double fUserMinTimeStep = -1
 

Detailed Description

Given a molecule G4DNAMoleculeEncounterStepper will calculate for its possible reactants what will be the minimum encounter time and the associated molecules.*

This model includes dynamical time steps as explained in "Computer-Aided Stochastic Modeling of the Radiolysis of Liquid Water", V. Michalik, M. Begusová, E. A. Bigildeev, Radiation Research, Vol. 149, No. 3 (Mar., 1998), pp. 224-236

Definition at line 61 of file G4DNAMoleculeEncounterStepper.hh.

Constructor & Destructor Documentation

G4DNAMoleculeEncounterStepper::G4DNAMoleculeEncounterStepper ( )

Definition at line 47 of file G4DNAMoleculeEncounterStepper.cc.

47  :
49  fMolecularReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable)),
50  fReactionModel(0)
51 {
52  fVerbose = 0;
53  fHasAlreadyReachedNullTime = false;
54 }
const G4ITReactionTable * fpReactionTable
G4DNAMoleculeEncounterStepper::~G4DNAMoleculeEncounterStepper ( )
virtual

Definition at line 66 of file G4DNAMoleculeEncounterStepper.cc.

67 {}
G4DNAMoleculeEncounterStepper::G4DNAMoleculeEncounterStepper ( const G4DNAMoleculeEncounterStepper right)

Definition at line 69 of file G4DNAMoleculeEncounterStepper.cc.

69  :
70  G4VITTimeStepper(right),
71  fMolecularReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
72 {
73  fVerbose = right.fVerbose ;
74  fMolecularReactionTable = right.fMolecularReactionTable;
75  fReactionModel = 0;
76  fHasAlreadyReachedNullTime = false;
77 }
const G4ITReactionTable * fpReactionTable

Member Function Documentation

G4double G4DNAMoleculeEncounterStepper::CalculateStep ( const G4Track trackA,
const G4double userMinTimeStep 
)
virtual

Implements G4VITTimeStepper.

Definition at line 87 of file G4DNAMoleculeEncounterStepper.cc.

References DBL_MAX, G4VITTimeStepper::fReactants, G4VITTimeStepper::fSampledMinTimeStep, G4VITTimeStepper::fUserMinTimeStep, G4BestUnit, G4cout, G4endl, GetMolecule(), G4Molecule::GetName(), and G4Track::GetTrackID().

88 {
89  // DEBUG
90  // G4cout << "G4MoleculeEncounterStepper::CalculateStep, time :" << G4ITTrackHolder::Instance()->GetGlobalTime() << G4endl;
91 
92  G4Molecule* moleculeA = GetMolecule(trackA);
93 
94 
95 #ifdef G4VERBOSE
96  if(fVerbose)
97  {
98  G4cout << "_______________________________________________________________________" << G4endl;
99  G4cout << "G4DNAMoleculeEncounterStepper::CalculateStep" << G4endl;
100  G4cout << "Incident molecule : " << moleculeA->GetName()
101  << " (" << trackA.GetTrackID() << ") "
102  << G4endl;
103  }
104 #endif
105 
106  //__________________________________________________________________
107  // Retrieve general informations for making reactions
108  const vector<const G4Molecule*>* reactivesVector =
109  fMolecularReactionTable -> CanReactWith(moleculeA);
110 
111  if(!reactivesVector)
112  {
113 #ifdef G4VERBOSE
114  // DEBUG
115  if(fVerbose > 1)
116  {
117  G4cout << "!!!!!!!!!!!!!!!!!!!!"<<G4endl;
118  G4cout << "!!! WARNING" << G4endl;
119  G4cout << "G4MoleculeEncounterStepper::CalculateStep will return infinity for the reaction because the molecule "
120  << moleculeA->GetName()
121  << " does not have any reactants given in the reaction table."
122  << G4endl;
123  G4cout << "!!!!!!!!!!!!!!!!!!!!"<<G4endl;
124  }
125 #endif
126  return DBL_MAX;
127  }
128 
129  G4int nbReactives = reactivesVector->size();
130 
131  if(nbReactives == 0)
132  {
133 #ifdef G4VERBOSE
134  // DEBUG
135  if(fVerbose)
136  {
137  G4cout << "!!!!!!!!!!!!!!!!!!!!"<<G4endl;
138  G4cout << "!!! WARNING" << G4endl;
139  G4cout << "G4MoleculeEncounterStepper::CalculateStep will return infinity for the reaction because the molecule "
140  << moleculeA->GetName()
141  << " does not have any reactants given in the reaction table."
142  << "This message can also result from a wrong implementation of the reaction table."
143  << G4endl;
144  G4cout << "!!!!!!!!!!!!!!!!!!!!"<<G4endl;
145  }
146 #endif
147  return DBL_MAX;
148  }
149  // DEBUG
150  // else
151  // {
152  // G4cout << "nb reactants : " << nbReactives << " pour mol "<< moleculeA -> GetName () << G4endl;
153  // for(int k=0 ; k < nbReactives ; k++)
154  // {
155  // G4cout << (*reactivesVector)[k]->GetName() << G4endl;
156  // }
157  // }
158 
159  fUserMinTimeStep = userMinTimeStep ;
160  if(fReactants) fReactants = 0 ;
161  fReactants = new vector<G4Track*>();
162 
164  fHasAlreadyReachedNullTime = false;
165 
166  fReactionModel -> Initialise(moleculeA, trackA) ;
167 
168  //__________________________________________________________________
169  // Start looping on possible reactants
170  for (G4int i=0 ; i<nbReactives ; i++)
171  {
172  const G4Molecule* moleculeB = (*reactivesVector)[i];
173 
174  //______________________________________________________________
175  // Retrieve reaction range
176  G4double R = -1 ; // reaction Range
177  R = fReactionModel -> GetReactionRadius(i);
178 
179  //______________________________________________________________
180  // Use KdTree algorithm to find closest reactants
182  -> FindNearest(moleculeA, moleculeB));
183 
184  RetrieveResults(trackA,moleculeA,moleculeB,R,results, true);
185  }
186 
187 #ifdef G4VERBOSE
188  // DEBUG
189  if(fVerbose)
190  {
191  G4cout << "G4MoleculeEncounterStepper::CalculateStep will finally return :"
192  << G4BestUnit(fSampledMinTimeStep, "Time")<< G4endl;
193 
194  if(fVerbose > 1)
195  {
196  G4cout << "TrackA: " << moleculeA->GetName() << " (" << trackA.GetTrackID() << ") can react with: ";
197 
198  vector<G4Track*>::iterator it;
199  for(it = fReactants->begin() ; it != fReactants->end() ; it++)
200  {
201  G4Track* trackB = *it;
202  G4cout << GetMolecule(trackB)->GetName() << " ("
203  << trackB->GetTrackID() << ") \t ";
204  }
205  G4cout << G4endl;
206  }
207  }
208 #endif
209  return fSampledMinTimeStep ;
210 }
static G4ThreadLocal G4double fUserMinTimeStep
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
G4TrackVectorHandle fReactants
const G4String & GetName() const
Definition: G4Molecule.cc:243
G4GLOB_DLL std::ostream G4cout
G4int GetTrackID() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
#define G4endl
Definition: G4ios.hh:61
G4double fSampledMinTimeStep
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4VDNAReactionModel * G4DNAMoleculeEncounterStepper::GetReactionModel ( )
inline

Definition at line 98 of file G4DNAMoleculeEncounterStepper.hh.

99 {
100  return fReactionModel;
101 }
void G4DNAMoleculeEncounterStepper::Prepare ( )
virtual

Reimplemented from G4VITTimeStepper.

Definition at line 79 of file G4DNAMoleculeEncounterStepper.cc.

References G4ITManager< T >::Instance(), and G4VITTimeStepper::Prepare().

80 {
81  // DEBUG
82  // G4cout << "G4DNAMoleculeEncounterStepper::PrepareForAllProcessors" << G4endl;
84  G4ITManager<G4Molecule>::Instance()->UpdatePositionMap();
85 }
static G4ITManager< T > * Instance()
virtual void Prepare()
void G4DNAMoleculeEncounterStepper::SetReactionModel ( G4VDNAReactionModel reactionModel)
inline

Definition at line 93 of file G4DNAMoleculeEncounterStepper.hh.

94 {
95  fReactionModel = reactionModel ;
96 }
void G4DNAMoleculeEncounterStepper::SetVerbose ( int  flag)
inline

Definition at line 103 of file G4DNAMoleculeEncounterStepper.hh.

104 {
105  fVerbose = flag;
106 }

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