Geant4  10.02.p03
G4DNASmoluchowskiDiffusion Class Reference

#include <G4DNASmoluchowskiDiffusion.hh>

Collaboration diagram for G4DNASmoluchowskiDiffusion:

Classes

struct  BoundingBox
 

Public Member Functions

 G4DNASmoluchowskiDiffusion (double epsilon=1e-5)
 
virtual ~G4DNASmoluchowskiDiffusion ()
 
double GetRandomDistance (double _time, double D)
 
double GetRandomTime (double distance, double D)
 
double EstimateCrossingTime (double proba, double distance, double D)
 
void PrepareReverseTable (double xmin, double xmax)
 
double GetInverseProbability (double proba)
 
double PlotInverse (double *x, double *)
 
double Plot (double *x, double *)
 
void InitialiseInverseProbability (double xmax=3e28)
 

Static Public Member Functions

static double ComputeS (double r, double D, double t)
 
static double ComputeDistance (double sTransform, double D, double t)
 
static double ComputeTime (double sTransform, double D, double r)
 
static double GetDifferential (double sTransform)
 
static double GetDensityProbability (double r, double _time, double D)
 
static double GetCumulativeProbability (double sTransform)
 

Public Attributes

std::vector< double > fInverse
 
int fNbins
 
double fEpsilon
 

Detailed Description

Definition at line 68 of file G4DNASmoluchowskiDiffusion.hh.

Constructor & Destructor Documentation

◆ G4DNASmoluchowskiDiffusion()

G4DNASmoluchowskiDiffusion::G4DNASmoluchowskiDiffusion ( double  epsilon = 1e-5)

Definition at line 50 of file G4DNASmoluchowskiDiffusion.cc.

50  : fEpsilon(epsilon)
51 {
52  fNbins = (int) trunc(1/fEpsilon);
53  // std::cout << "fNbins: " << fNbins << std::endl;
54 #ifdef DNADEV
55  assert(fNbins > 0);
56 #endif
57  fInverse.resize(fNbins+2); // trunc sous-estime + borne max a rajouter ==> 2
58 
59  // std::cout << "fInverse.capacity(): "<< fInverse.capacity() << std::endl;
60 }
double epsilon(double density, double temperature)
Here is the caller graph for this function:

◆ ~G4DNASmoluchowskiDiffusion()

G4DNASmoluchowskiDiffusion::~G4DNASmoluchowskiDiffusion ( )
virtual

Definition at line 62 of file G4DNASmoluchowskiDiffusion.cc.

63 {
64 }
Here is the call graph for this function:

Member Function Documentation

◆ ComputeDistance()

static double G4DNASmoluchowskiDiffusion::ComputeDistance ( double  sTransform,
double  D,
double  t 
)
inlinestatic

Definition at line 80 of file G4DNASmoluchowskiDiffusion.hh.

81  {
82  return sTransform * 2. * std::sqrt(D * t);
83  }
double D(double temp)
Here is the caller graph for this function:

◆ ComputeS()

static double G4DNASmoluchowskiDiffusion::ComputeS ( double  r,
double  D,
double  t 
)
inlinestatic

Definition at line 74 of file G4DNASmoluchowskiDiffusion.hh.

75  {
76  double sTransform = r / (2. * std::sqrt(D * t));
77  return sTransform;
78  }
double D(double temp)

◆ ComputeTime()

static double G4DNASmoluchowskiDiffusion::ComputeTime ( double  sTransform,
double  D,
double  r 
)
inlinestatic

Definition at line 85 of file G4DNASmoluchowskiDiffusion.hh.

86  {
87  return std::pow(r / sTransform, 2.) / (4. * D);
88  }
double D(double temp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EstimateCrossingTime()

double G4DNASmoluchowskiDiffusion::EstimateCrossingTime ( double  proba,
double  distance,
double  D 
)
inline

Definition at line 108 of file G4DNASmoluchowskiDiffusion.hh.

111  {
112  double sTransform = GetInverseProbability(proba);
113  return ComputeTime(sTransform, D, distance);
114  }
static double ComputeTime(double sTransform, double D, double r)
double D(double temp)
Here is the call graph for this function:

◆ GetCumulativeProbability()

static double G4DNASmoluchowskiDiffusion::GetCumulativeProbability ( double  sTransform)
inlinestatic

Definition at line 302 of file G4DNASmoluchowskiDiffusion.hh.

303  {
304  static double constant = 2./std::sqrt(3.141592653589793);
305  return erfc(sTransform) + constant*sTransform*std::exp(-sTransform*sTransform);
306  }
Here is the caller graph for this function:

◆ GetDensityProbability()

static double G4DNASmoluchowskiDiffusion::GetDensityProbability ( double  r,
double  _time,
double  D 
)
inlinestatic

Definition at line 127 of file G4DNASmoluchowskiDiffusion.hh.

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  }
double D(double temp)
Here is the caller graph for this function:

◆ GetDifferential()

static double G4DNASmoluchowskiDiffusion::GetDifferential ( double  sTransform)
inlinestatic

Definition at line 121 of file G4DNASmoluchowskiDiffusion.hh.

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  }
Here is the caller graph for this function:

◆ GetInverseProbability()

double G4DNASmoluchowskiDiffusion::GetInverseProbability ( double  proba)
inline

Definition at line 308 of file G4DNASmoluchowskiDiffusion.hh.

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  }
Here is the caller graph for this function:

◆ GetRandomDistance()

double G4DNASmoluchowskiDiffusion::GetRandomDistance ( double  _time,
double  D 
)
inline

Definition at line 92 of file G4DNASmoluchowskiDiffusion.hh.

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  }
static double ComputeDistance(double sTransform, double D, double t)
#define G4UniformRand()
Definition: Randomize.hh:97
double D(double temp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetRandomTime()

double G4DNASmoluchowskiDiffusion::GetRandomTime ( double  distance,
double  D 
)
inline

Definition at line 101 of file G4DNASmoluchowskiDiffusion.hh.

102  {
103  double proba = G4UniformRand();
104  double sTransform = GetInverseProbability(proba);
105  return ComputeTime(sTransform, D, distance);
106  }
#define G4UniformRand()
Definition: Randomize.hh:97
static double ComputeTime(double sTransform, double D, double r)
double D(double temp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseInverseProbability()

void G4DNASmoluchowskiDiffusion::InitialiseInverseProbability ( double  xmax = 3e28)
inline

Definition at line 355 of file G4DNASmoluchowskiDiffusion.hh.

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  }
void PrepareReverseTable(double xmin, double xmax)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Plot()

double G4DNASmoluchowskiDiffusion::Plot ( double *  x,
double *   
)
inline

Definition at line 349 of file G4DNASmoluchowskiDiffusion.hh.

350  {
351  return GetDifferential(x[0]);
352  }
static double GetDifferential(double sTransform)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PlotInverse()

double G4DNASmoluchowskiDiffusion::PlotInverse ( double *  x,
double *   
)
inline

Definition at line 344 of file G4DNASmoluchowskiDiffusion.hh.

345  {
346  return GetInverseProbability(x[0]);
347  }
Here is the call graph for this function:

◆ PrepareReverseTable()

void G4DNASmoluchowskiDiffusion::PrepareReverseTable ( double  xmin,
double  xmax 
)
inline

Definition at line 261 of file G4DNASmoluchowskiDiffusion.hh.

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  }
Int_t index
static double GetCumulativeProbability(double sTransform)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fEpsilon

double G4DNASmoluchowskiDiffusion::fEpsilon

Definition at line 370 of file G4DNASmoluchowskiDiffusion.hh.

◆ fInverse

std::vector<double> G4DNASmoluchowskiDiffusion::fInverse

Definition at line 368 of file G4DNASmoluchowskiDiffusion.hh.

◆ fNbins

int G4DNASmoluchowskiDiffusion::fNbins

Definition at line 369 of file G4DNASmoluchowskiDiffusion.hh.


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