Geant4  10.02
G4DNASmoluchowskiDiffusion.hh
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 /*
27  * G4DNASmoluchowskiDiffusion.hh
28  *
29  * Created on: 2 févr. 2015
30  * Author: matkara
31  */
32 
33 #ifndef G4DNASMOLUCHOWSKIDIFFUSION_HH_
34 #define G4DNASMOLUCHOWSKIDIFFUSION_HH_
35 
36 //#if __cplusplus >= 201103L
37 #include <cstdlib>
38 #include <cmath>
39 #include <vector>
40 #include <iostream>
41 #include <algorithm>
42 
43 //#define DNADEV_TEST
44 
45 #ifdef DNADEV_TEST
46 #include <TF1.h>
47 #endif
48 
49 #include <cassert>
50 
51 #ifndef DNADEV_TEST
52 #include "globals.hh"
53 #include "Randomize.hh"
54 #endif
55 
56 #ifdef DNADEV_TEST
57 #include "TRandom.h"
58 TRandom root_random;
59 double G4UniformRand()
60 {
61  return root_random.Rndm();
62 }
63 
64 #define G4cout std::cout
65 #define G4endl std::endl
66 #endif
67 
69 {
70 public:
71  G4DNASmoluchowskiDiffusion(double epsilon = 1e-5);
73 
74  static double ComputeS(double r, double D, double t)
75  {
76  double sTransform = r / (2. * std::sqrt(D * t));
77  return sTransform;
78  }
79 
80  static double ComputeDistance(double sTransform, double D, double t)
81  {
82  return sTransform * 2. * std::sqrt(D * t);
83  }
84 
85  static double ComputeTime(double sTransform, double D, double r)
86  {
87  return std::pow(r / sTransform, 2.) / (4. * D);
88  }
89 
90  //====================================================
91 
92  double GetRandomDistance(double _time, double D)
93  {
94  double proba = G4UniformRand();
95 // G4cout << "proba = " << proba << G4endl;
96  double sTransform = GetInverseProbability(proba);
97 // G4cout << "sTransform = " << sTransform << G4endl;
98  return ComputeDistance(sTransform, D, _time);
99  }
100 
101  double GetRandomTime(double distance, double D)
102  {
103  double proba = G4UniformRand();
104  double sTransform = GetInverseProbability(proba);
105  return ComputeTime(sTransform, D, distance);
106  }
107 
108  double EstimateCrossingTime(double proba,
109  double distance,
110  double D)
111  {
112  double sTransform = GetInverseProbability(proba);
113  return ComputeTime(sTransform, D, distance);
114  }
115 
116  //====================================================
117  // 1-value transformation
118 
119  // WARNING : this is NOT the differential probability
120  // this is the derivative of the function GetCumulativeProbability
121  static double GetDifferential(double sTransform)
122  {
123  static double constant = -4./std::sqrt(3.141592653589793);
124  return sTransform*sTransform*std::exp(-sTransform*sTransform)*constant; // -4*sTransform*sTransform*exp(-sTransform*sTransform)/sqrt(3.141592653589793)
125  }
126 
127  static double GetDensityProbability(double r, double _time, double D)
128  {
129  static double my_pi = 3.141592653589793;
130  static double constant = 4.*my_pi/std::pow(4.*my_pi, 1.5);
131  return r*r/std::pow(D * _time,1.5)*std::exp(-r*r/(4. * D * _time))*constant;
132  }
133 
134  //====================================================
135  // BOUNDING BOX
136  struct BoundingBox
137  {
138  double fXmax;
139  double fXmin;
140  double fXmaxDef;
141  double fXminDef;
142  double fToleranceY;
143  double fSum;
145 
147  {
151  };
152 
154 
155  BoundingBox(double xmin,
156  double xmax,
157  double toleranceY) :
158  fXmax(xmax), fXmin(xmin),
159  fToleranceY(toleranceY),
160  fSum(0)
161  {
162  if(fXmax < fXmin)
163  {
164  double tmp = fXmin;
165  fXmin = fXmax;
166  fXmax = tmp;
167  }
168 
169  fXminDef = fXmin;
170  fXmaxDef = fXmax;
171  fPreviousAction = BoundingBox::Undefined;
172  fIncreasingCumulativeFunction = (GetCumulativeProbability(fXmax) - GetCumulativeProbability(fXmin))/(fXmax-fXmin);
173  }
174 
175  void Print()
176  {
177  G4cout << "fXmin: " << fXmin << " | fXmax: " << fXmax << G4endl;
178  }
179 
180  bool Propose(double proposedXValue,
181  double proposedProba,
182  double nextProba,
183  double& returnedValue)
184  {
185 // G4cout << "---------------------------" << G4endl;
186 // G4cout << "Proposed x value: " << proposedXValue
187 // << "| proposedProba: " << proposedProba
188 // << "| nextProba: " << nextProba
189 // << " | fXmin: " << fXmin << " (" << G4DNASmoluchowskiDiffusion::GetCumulativeProbability(fXmin) <<")"
190 // << " | fXmax: " << fXmax << " (" << G4DNASmoluchowskiDiffusion::GetCumulativeProbability(fXmax) <<")"
191 // << G4endl;
192 
193  bool returnFlag = false;
194 
195  if(proposedProba < nextProba-fToleranceY) // proba trop petite ==> augmente
196  {
197  // G4cout << "proposedProba < nextProba-fToleranceY" << G4endl;
198 
199  if(fIncreasingCumulativeFunction > 0) // croissant
200  {
201  if(proposedXValue > fXmin)
202  fXmin = proposedXValue;
203  }
204  else if(fIncreasingCumulativeFunction < 0) // decroissant
205  {
206  if(proposedXValue < fXmax)
207  fXmax = proposedXValue;
208  }
209 
210  returnedValue = (fXmax + fXmin)/2;
211  returnFlag = false;
212  fPreviousAction = BoundingBox::IncreaseProba;
213  }
214  else if(proposedProba > nextProba+fToleranceY) // proba trop grande
215  {
216  // G4cout << "proposedProba > nextProba+fToleranceY" << G4endl;
217 
218  if(fIncreasingCumulativeFunction>0)
219  {
220  if(proposedXValue < fXmax)
221  fXmax = proposedXValue;
222  }
223  else if(fIncreasingCumulativeFunction<0)
224  {
225  if(proposedXValue > fXmin)
226  {
227  fXmin = proposedXValue;
228  }
229  }
230 
231  returnedValue = (fXmax + fXmin)/2;
232  returnFlag = false;
233  fPreviousAction = BoundingBox::DecreaseProba;
234  }
235  else
236  {
237  // G4cout << "IN THE INTERVAL !! : " << nextProba << G4endl;
238  fSum = proposedProba;
239 
240  // Assuming search for next proba is increasing
241  if(fIncreasingCumulativeFunction<0)
242  {
243  fXmin = fXminDef;
244  fXmax = proposedXValue;
245  }
246  else if(fIncreasingCumulativeFunction>0)
247  {
248  fXmin = proposedXValue;
249  fXmax = fXmaxDef;
250  }
251  returnFlag = true;
252  fPreviousAction = BoundingBox::Undefined;
253  }
254 
255  return returnFlag;
256  }
257  };
258  // END OF BOUNDING BOX
259  //==============================
260 
261  void PrepareReverseTable(double xmin, double xmax)
262  {
263  double x = xmax;
264  int index = 0;
265  double nextProba = fEpsilon;
266  double proposedX;
267 
268  BoundingBox boundingBox(xmin, xmax, fEpsilon*1e-5);
269 
270  while(index <= fNbins)
271  // in case GetCumulativeProbability is exact (digitally speaking), replace with:
272  // while(index <= fNbins+1)
273  {
274  nextProba = fEpsilon*index;
275 
276  double newProba = GetCumulativeProbability(x);
277 
278  if(boundingBox.Propose(x, newProba, nextProba, proposedX))
279  {
280  fInverse[index] = x;
281  index++;
282  }
283  else
284  {
285  if(x == proposedX)
286  {
287  G4cout << "BREAK : x= " << x << G4endl;
288  G4cout << "index= " << index << G4endl;
289  G4cout << "nextProba= " << nextProba << G4endl;
290  G4cout << "newProba= " << newProba << G4endl;
291  abort();
292  }
293  x = proposedX;
294  }
295  }
296 
297  fInverse[fNbins+1] = 0; // P(1) = 0, because we want it exact !
298  // Tips to improve the exactness: get an better value of pi, get better approximation of erf and exp, use long double instead of double
299 // boundingBox.Print();
300  }
301 
302  static double GetCumulativeProbability(double sTransform)
303  {
304  static double constant = 2./std::sqrt(3.141592653589793);
305  return erfc(sTransform) + constant*sTransform*std::exp(-sTransform*sTransform);
306  }
307 
308  double GetInverseProbability(double proba) // returns sTransform
309  {
310  size_t index_low = (size_t) trunc(proba/fEpsilon);
311 
312  if(index_low == 0) // assymptote en 0
313  {
314  index_low = 1;
315  size_t index_up = 2;
316  double low_y = fInverse[index_low];
317  double up_y = fInverse[index_up];
318  double low_x = index_low*fEpsilon;
319  double up_x = proba+fEpsilon;
320  double tangente = (low_y-up_y)/(low_x - up_x); // ou utiliser GetDifferential(proba) ?
321  // double tangente = GetDifferential(proba);
322  return low_y + tangente*(proba-low_x);
323  }
324 
325  size_t index_up = index_low+1;
326  if(index_low > fInverse.size()) return fInverse.back();
327  double low_y = fInverse[index_low];
328  double up_y = fInverse[index_up];
329 
330  double low_x = index_low*fEpsilon;
331  double up_x = low_x+fEpsilon;
332 
333  if(up_x > 1) // P(1) = 0
334  {
335  up_x = 1;
336  up_y = 0; // more general : fInverse.back()
337  }
338 
339  double tangente = (low_y-up_y)/(low_x - up_x);
340 
341  return low_y + tangente*(proba-low_x);
342  }
343 
344  double PlotInverse(double* x, double* )
345  {
346  return GetInverseProbability(x[0]);
347  }
348 
349  double Plot(double* x, double* )
350  {
351  return GetDifferential(x[0]);
352  }
353 
354 
355  void InitialiseInverseProbability(double xmax = 3e28)
356  {
357  // x > x'
358  // P'(x) = p(x') = lim(x->x') (P(x') - P(x))/(x'-x)
359  // x'-x = (P(x') - P(x))/p(x')
360  // x = x' - (P(x') - P(x))/p(x')
361 
362  // fInverse initialized in the constructor
363 
364  assert(fNbins !=0);
365  PrepareReverseTable(0,xmax);
366  }
367 
368  std::vector<double> fInverse;
369  int fNbins;
370  double fEpsilon;
371 };
372 
373 #endif /* SOURCE_PROCESSES_ELECTROMAGNETIC_DNA_MODELS_G4DNASMOLUCHOWSKIDIFFUSION_HH_ */
void InitialiseInverseProbability(double xmax=3e28)
static double ComputeDistance(double sTransform, double D, double t)
double EstimateCrossingTime(double proba, double distance, double D)
double PlotInverse(double *x, double *)
static double GetCumulativeProbability(double sTransform)
void PrepareReverseTable(double xmin, double xmax)
G4DNASmoluchowskiDiffusion(double epsilon=1e-5)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
BoundingBox(double xmin, double xmax, double toleranceY)
double GetRandomTime(double distance, double D)
static double ComputeTime(double sTransform, double D, double r)
static double ComputeS(double r, double D, double t)
const G4double x[NPOINTSGL]
static double GetDifferential(double sTransform)
double D(double temp)
#define G4endl
Definition: G4ios.hh:61
static double GetDensityProbability(double r, double _time, double D)
bool Propose(double proposedXValue, double proposedProba, double nextProba, double &returnedValue)
double epsilon(double density, double temperature)
double Plot(double *x, double *)
double GetRandomDistance(double _time, double D)