Geant4  10.00.p02
G4DNABrownianTransportation.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4DNABrownianTransportation.cc 82326 2014-06-16 09:19:18Z gcosmo $
27 //
28 // Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29 //
30 // WARNING : This class is released as a prototype.
31 // It might strongly evolve or even disapear in the next releases.
32 //
33 // History:
34 // -----------
35 // 10 Oct 2011 M.Karamitros created
36 //
37 // -------------------------------------------------------------------
38 
41 
42 #include <CLHEP/Random/Stat.h>
43 
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "Randomize.hh"
48 #include "G4Molecule.hh"
49 #include "G4RandomDirection.hh"
50 #include "G4ParticleTable.hh"
51 #include "G4SafetyHelper.hh"
53 #include "G4UnitsTable.hh"
54 #include "G4NistManager.hh"
56 
57 using namespace std;
58 
59 #ifndef State
60 #define State(theXInfo) (GetState<G4ITBrownianState>()->theXInfo)
61 #endif
62 
63 //#ifndef State
64 //#define State(theXInfo) (fTransportationState->theXInfo)
65 //#endif
66 
67 //COLOR FOR DEBUGING
68 //#define RED_ON_WHITE "\033[0;31m"
69 //#define GREEN "\033[32;40m"
70 #define GREEN_ON_BLUE "\033[1;32;44m"
71 #define RESET "\033[0m"
72 
73 static double InvErf(double x)
74 {
75  return CLHEP::HepStat::inverseErf(x);
76 }
77 
78 static double InvErfc(double x)
79 {
80  return CLHEP::HepStat::inverseErf(1.-x);
81 }
82 
84  G4ITTransportation(aName, verbosity)
85  {
86  //ctor
87  fpState.reset(new G4ITBrownianState());
89  verboseLevel = 1;
92  fpWaterDensity = 0;
93 }
94 
96 {;}
97 
99  G4ITTransportation(right)
100 {
101  //copy ctor
102  SetProcessSubType(61);
103  verboseLevel = right.verboseLevel;
105  fNistWater = right.fNistWater;
107 }
108 
110 {
111  if (this == &rhs) return *this; // handle self assignment
112  //assignment operator
113  return *this;
114 }
115 
117 {
118  fPathLengthWasCorrected = false;
119 }
120 
122 {
123  fpState.reset(new G4ITBrownianState());
126 }
127 
129 {
130  if(verboseLevel > 0)
131  {
132  G4cout << G4endl << GetProcessName() << ": for "
133  << setw(24) << particle.GetParticleName()
134  << "\tSubType= " << GetProcessSubType() << G4endl;
135  }
136  // Initialize water density pointer
138 }
139 
141  const G4Step& step,
142  const double timeStep,
143  double& spaceStep)
144 {
145  // G4cout << "G4ITBrownianTransportation::ComputeStep" << G4endl;
146 
147  /* If this method is called, this step
148  * cannot be geometry limited.
149  * In case the step is limited by the geometry,
150  * this method should not be called.
151  * The fTransportEndPosition calculated in
152  * the method AlongStepIL should be taken
153  * into account.
154  * In order to do so, the flag IsLeadingStep
155  * is on. Meaning : this track has the minimum
156  * interaction length over all others.
157  */
158  if(GetIT(track)->GetTrackingInfo()->IsLeadingStep())
159  {
160  const G4VITProcess* ITProc = ((const G4VITProcess*) step.GetPostStepPoint()->GetProcessDefinedStep());
161  bool makeException = true;
162 
163  if(ITProc && ITProc->ProposesTimeStep()) makeException = false;
164 
165  if(makeException)
166  {
167 
168  G4ExceptionDescription exceptionDescription ;
169  exceptionDescription << "ComputeStep is called while the track has the minimum interaction time";
170  exceptionDescription << " so it should not recompute a timeStep ";
171  G4Exception("G4DNABrownianTransportation::ComputeStep","G4DNABrownianTransportation001",
172  FatalErrorInArgument,exceptionDescription);
173  }
174  }
175 
176  State(fGeometryLimitedStep) = false;
177  // TODO : generalize this process to all kind of brownian objects
178  // G4ITBrownianObject* ITBrown = GetITBrownianObject(track) ;
179  // G4double diffCoeff = ITBrown->GetDiffusionCoefficient(track.GetMaterial());
180  G4Molecule* molecule = GetMolecule(track) ;
181  G4double diffCoeff = molecule->GetDiffusionCoefficient();
182 
183  if(timeStep > 0)
184  {
185  spaceStep = DBL_MAX;
186 
187  while(spaceStep > State(endpointDistance))
188  // Probably inefficient when the track is close close to boundaries
189  // it goes with fUserMaximumTimeBeforeReachingBoundary == false
190  // fUserMaximumTimeBeforeReachingBoundary == true, it should never loop
191  {
192  G4double x = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
193  G4double y = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
194  G4double z = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
195 
196  spaceStep = sqrt(x*x + y*y + z*z);
197  }
198  // State(fTransportEndPosition).set(x + track.GetPosition().x(),
199  // y + track.GetPosition().y(),
200  // z + track.GetPosition().z());
201 
202  State(fTransportEndPosition)= spaceStep*step.GetPostStepPoint()->GetMomentumDirection() + track.GetPosition();
203  }
204  else
205  {
206  spaceStep = 0. ;
207  State(fTransportEndPosition) = track.GetPosition() ;
208  }
209 
210  State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime() + timeStep ;
211  State(fEndGlobalTimeComputed) = true ;
212 
213 #ifdef G4VERBOSE
214  // DEBUG
215  if(fVerboseLevel>1)
216  {
218  << "G4ITBrownianTransportation::ComputeStep() : "
219  << " trackID : " << track.GetTrackID()
220  << " : Molecule name: " << molecule-> GetName()
221  << G4endl
222  << "Diffusion length : " << G4BestUnit(spaceStep, "Length")
223  << " within time step : " << G4BestUnit(timeStep,"Time")
224  << RESET
225  << G4endl<< G4endl;
226  }
227 #endif
228 }
229 
231 {
233 
234 #ifdef G4VERBOSE
235  // DEBUG
236  if(fVerboseLevel>1)
237  {
239  << "G4ITBrownianTransportation::PostStepDoIt() :"
240  << " trackID : " << track.GetTrackID()
241  << " Molecule name: " << GetMolecule(track)-> GetName() << G4endl;
242  G4cout<< "Diffusion length : "<< G4BestUnit(step.GetStepLength(),"Length") <<" within time step : "
243  << G4BestUnit(step.GetDeltaTime(),"Time") << "\t"
244  << " Current global time : " << G4BestUnit(track.GetGlobalTime(),"Time")
245  << RESET
246  << G4endl<< G4endl;
247  }
248 #endif
249  return &fParticleChange ;
250 }
251 
253  const G4Track& track)
254 {
255 
256 #ifdef G4VERBOSE
257  // DEBUG
258  if (fVerboseLevel>1)
259  {
261  << setw(18)<< "G4DNABrownianTransportation::Diffusion :"
262  << setw(8) << GetIT(track)->GetName()
263  << "\t trackID:" << track.GetTrackID() <<"\t"
264  << " Global Time = " << G4BestUnit(track.GetGlobalTime(),"Time")
265  << RESET
266  << G4endl<< G4endl;
267  }
268 #endif
269 
270  G4Material* material = track.GetMaterial();
271 // if (material != fNistWater && material->GetBaseMaterial() != fNistWater)
272 
273  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
274 
275  if(waterDensity == 0.0)
276 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
277  {
278  G4cout << "A track is outside water material : trackID"<< track.GetTrackID() << " (" << GetMolecule(track)->GetName() <<")" << G4endl;
279  G4cout << "Local Time : "<< (track.GetLocalTime()) /s<<G4endl;
280  G4cout<< "Step Number :" << track.GetCurrentStepNumber() <<G4endl;
281 
284 
285  // Got pb with :
286  // fParticleChange.ProposeTrackStatus(fStopAndKill);
287  // It makes the tracks going straight without killing them
288 
289  return ; // &fParticleChange is the final returned object
290  }
291 
292  G4double costheta = (2*G4UniformRand()-1);
293  G4double theta = acos (costheta);
294  G4double phi = 2*pi*G4UniformRand();
295 
296  G4double xMomentum = cos(phi)* sin(theta);
297  G4double yMomentum = sin(theta)*sin(phi);
298  G4double zMomentum = costheta;
299 
300  fParticleChange.ProposeMomentumDirection(xMomentum, yMomentum, zMomentum);
301  State(fMomentumChanged) = true;
303 
304  // G4cout << "BROWN : Propose new direction :" << G4ThreeVector(xMomentum, yMomentum, zMomentum) << G4endl;
305 
306  // Alternative
307  //fParticleChange.ProposeMomentumDirection(G4RandomDirection());
308 
309  return; // &fParticleChange is the final returned object
310 }
311 
312 
314  const G4Track& track,
315  G4double previousStepSize,
316  G4double currentMinimumStep,
317  G4double& currentSafety,
318  G4GPILSelection* selection)
319 {
321  previousStepSize,
322  currentMinimumStep,
323  currentSafety,
324  selection);
325 
326  G4double diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
327  // G4double diffusionCoefficient = GetITBrownianObject(track)->GetDiffusionCoefficient(track.GetMaterial());
328 
329  if(State(fGeometryLimitedStep))
330  {
331  // 99 % of the space step distribution is lower than
332  // d_99 = 8 * sqrt(D*t)
333  // where t is the corresponding time step
334  // so by inversion :
336  {
337  State(theInteractionTimeLeft) = (geometryStepLength*geometryStepLength)/(64 * diffusionCoefficient);
338  }
339  else // Will use a random time
340  {
341  State(theInteractionTimeLeft) = 1/(4*diffusionCoefficient) * pow(geometryStepLength/InvErfc(G4UniformRand()),2);
342  }
343 
344  State(fCandidateEndGlobalTime) = track.GetGlobalTime() + State(theInteractionTimeLeft);
345  State(fPathLengthWasCorrected) = false;
346  }
347  else
348  {
349  geometryStepLength = 2*sqrt(diffusionCoefficient*State(theInteractionTimeLeft) ) *InvErf(G4UniformRand());
350  State(fPathLengthWasCorrected) = true;
351  State(endpointDistance) = geometryStepLength;
352  }
353 
354  return geometryStepLength ;
355 }
356 
358 //
359 // Initialize ParticleChange (by setting all its members equal
360 // to corresponding members in G4Track)
362  const G4Step& step )
363 {
365  Diffusion(track);
366  return &fParticleChange;
367 }
ThreeVector shoot(const G4int Ap, const G4int Af)
const std::vector< double > * GetDensityTableFor(const G4Material *) const
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
G4double GetLocalTime() const
G4ParticleChangeForTransport fParticleChange
G4DNABrownianTransportation & operator=(const G4DNABrownianTransportation &other)
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetStepLength() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
G4double z
Definition: TRTMaterials.hh:39
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &)
size_t GetIndex() const
Definition: G4Material.hh:260
const G4double pi
const G4ThreeVector & GetPosition() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
G4double GetDiffusionCoefficient() const
Returns the diffusion coefficient D.
Definition: G4Molecule.cc:389
virtual void StartTracking(G4Track *aTrack)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual const G4String & GetName() const =0
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
int G4int
Definition: G4Types.hh:78
static double InvErfc(double x)
static G4NistManager * Instance()
const G4String & GetParticleName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
static const double s
Definition: G4SIunits.hh:150
G4StepPoint * GetPreStepPoint() const
#define State(theXInfo)
const G4String & GetName() const
Returns the name of the molecule.
Definition: G4Molecule.cc:243
virtual void StartTracking(G4Track *aTrack)
const G4ThreeVector & GetMomentumDirection() const
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
void SetInstantiateProcessState(G4bool flag)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4double GetDeltaTime() const
G4int GetTrackID() const
G4double GetGlobalTime() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
virtual void ComputeStep(const G4Track &, const G4Step &, const double, double &)
void ProposeEnergy(G4double finalEnergy)
G4::shared_ptr< G4ProcessState > fpState
void Diffusion(const G4Track &track)
G4DNABrownianTransportation(const G4String &aName="DNABrownianTransportation", G4int verbosityLevel=1)
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
#define GREEN_ON_BLUE
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
Class Description The dynamic molecule holds all the data that change for a molecule It has a pointer...
Definition: G4Molecule.hh:76
static double InvErf(double x)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4bool ProposesTimeStep() const
G4VITProcess inherits from G4VProcess.
Definition: G4VITProcess.hh:74
void SetMomentumChanged(G4bool b)
#define DBL_MAX
Definition: templates.hh:83
G4double * theInteractionTimeLeft
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
const std::vector< G4double > * fpWaterDensity
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
G4GPILSelection