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