ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNABrownianTransportation.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4DNABrownianTransportation.hh
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 //
27 // Author: Mathieu Karamitros
28 
29 // The code is developed in the framework of the ESA AO7146
30 //
31 // We would be very happy hearing from you, send us your feedback! :)
32 //
33 // In order for Geant4-DNA to be maintained and still open-source,
34 // article citations are crucial.
35 // If you use Geant4-DNA chemistry and you publish papers about your software,
36 // in addition to the general paper on Geant4-DNA:
37 //
38 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39 //
40 // we would be very happy if you could please also cite the following
41 // reference papers on chemistry:
42 //
43 // J. Comput. Phys. 274 (2014) 841-882
44 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45 
46 #ifndef G4ITBROWNIANTRANSPORTATION_H
47 #define G4ITBROWNIANTRANSPORTATION_H
48 
49 #include "G4ITTransportation.hh"
50 
51 class G4SafetyHelper;
52 class G4Molecule;
53 
54 // experimental
56 {
57 public:
59  virtual ~G4BrownianAction(){;}
60 
61 // virtual G4double GetDiffusionCoefficient(G4Material*,
62 // G4Molecule*) { return 0;}
63 
64  // If returns true: track is killed
65  virtual void Transport(const G4Track&,
67 };
68 
69 
70 /* \brief {The transportation method implemented is the one from
71  * Ermak-McCammon : J. Chem. Phys. 69, 1352 (1978).
72  * To compute time and space intervals to reach a volume boundary,
73  * there are two alternative methods proposed by this process.
74  *
75  * ** Method 1 selects a minimum distance to the next
76  * boundary using to the following formula:
77  *
78  * --> t_min = (safety* safety) / (8 * diffusionCoefficient);
79  * this corresponds to 5% probability of the Brownian particle to cross
80  * the boundary - isotropic distance to nearest boundary (safety) is used
81  *
82  * OR if the flag "speed me up" is on:
83  *
84  * --> t_min = (geometryStepLength * geometryStepLength) / diffusionCoefficient;
85  * this corresponds to 50% probability of the Brownian particle to cross
86  * the boundary - distance along current direction to nearest boundary is used
87  *
88  * NB: By default, method 1 with the flag "speed me up is used".
89  * In addition, one may want to used the minimum time step limit defined
90  * in G4Scheduler through the G4UserTimeStepAction. If so, speed level might
91  * be set to 2. But minimum time steps have to be set in the user class.
92  *
93  * ** Method 2 can randomly compute the time to the next boundary using the
94  * following formula:
95  *
96  * t_random = 1 / (4 * diffusionCoefficient)* pow(geometryStepLength /
97  * InvErfc(G4UniformRand()),2);
98  * For release 10.1, this is using the 1D cumulative density function.
99  *
100  * At each diffusion step, the direction of the particle is randomly selected.
101  * For now, the geometryStepLength corresponds to the distance to the
102  * nearest boundary along the direction of diffusion which selected randomly.
103  *
104  * Method 2 is currently deactivated by default.
105  * }
106  */
107 
109 {
110 public:
112  "DNABrownianTransportation",
113  G4int verbosityLevel = 0);
118 
119  inline void SetBrownianAction(G4BrownianAction*);
120 
121  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
122 
123  virtual void StartTracking(G4Track* aTrack);
124 
125  virtual void ComputeStep(const G4Track&,
126  const G4Step&,
127  const double,
128  double&);
129 
130  virtual G4double
132  G4double /*previousStepSize*/,
133  G4double /*currentMinimumStep*/,
134  G4double& /*currentSafety*/,
135  G4GPILSelection* /*selection*/);
136 
137  virtual G4VParticleChange* PostStepDoIt(const G4Track& track, const G4Step&);
138 
139  virtual G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step&);
140 
141  // Boundary is crossed at time at which:
142  // * either 5% of the distribution might be over boundary - the last position
143  // is adjusted on boundary
144  // * or if speedUp (from level 1) is activated - 50% of the distribution might
145  // be over boundary, the particles are also allowed to jump over boundary
146  inline void UseMaximumTimeBeforeReachingBoundary(bool flag = true)
147  {
149  }
150 
151  // Random sampling time at which boundary is crossed
152  // WARNING: For release 10.1, this is a 1D approximation for sampling time
153  // but 3D for diffusion
154  // If speed up IS activated, particles are allowed jump over barrier
155  inline void UseCumulativeDensitFunction(bool flag = true)
156  {
157  if(flag == true)
158  {
160  return;
161  }
163  }
164 
165  // Use limiting time steps defined in the scheduler
166  inline void UseLimitingTimeSteps(bool flag = true)
167  {
169  }
170 
171  inline void SpeedLevel(int level)
172  {
173  if(level < 0) level =0;
174  else if(level > 2) level = 2;
175 
176  switch(level)
177  {
178  case 0:
179  fSpeedMeUp = false;
181  return;
182 
183  case 1:
184  fSpeedMeUp = true;
186  return;
187 
188  case 2:
189  //======================================================================
190  // NB: BE AWARE THAT IF NO MIN TIME STEPS NO TIME STEPS HAVE BEEN
191  // PROVIDED TO G4Scheduler THIS LEVEL MIGHT BE SLOWER THAN LEVEL 1
192  //======================================================================
193  fSpeedMeUp = true;
195  return;
196  }
197  }
198 
199 protected:
200 
201  G4double ComputeGeomLimit(const G4Track& track,
202  G4double& presafety,
203  G4double limit);
204 
205  void Diffusion(const G4Track& track);
206 
207  //________________________________________________________________
208  // Process information
209  struct G4ITBrownianState : public G4ITTransportationState
210  {
211  public:
214  {
215  ;
216  }
217  virtual G4String GetType()
218  {
219  return "G4ITBrownianState";
220  }
221 
226  };
227 
230 
234 
235  // Water density table
236  const std::vector<G4double>* fpWaterDensity;
237 
239 };
240 
241 
243 {
244  fpBrownianAction = brownianAction;
245 }
246 
247 #endif // G4ITBROWNIANTRANSPORTATION_H