Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNASmoluchowskiDiffusion Class Reference

#include <G4DNASmoluchowskiDiffusion.hh>

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 70 of file G4DNASmoluchowskiDiffusion.hh.

Constructor & Destructor Documentation

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 }
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
double epsilon(double density, double temperature)

Here is the call graph for this function:

G4DNASmoluchowskiDiffusion::~G4DNASmoluchowskiDiffusion ( )
virtual

Definition at line 62 of file G4DNASmoluchowskiDiffusion.cc.

63 {
64 }

Member Function Documentation

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

Definition at line 82 of file G4DNASmoluchowskiDiffusion.hh.

83  {
84  return sTransform * 2. * std::sqrt(D * t);
85  }
double D(double temp)

Here is the caller graph for this function:

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

Definition at line 76 of file G4DNASmoluchowskiDiffusion.hh.

77  {
78  double sTransform = r / (2. * std::sqrt(D * t));
79  return sTransform;
80  }
double D(double temp)
static double G4DNASmoluchowskiDiffusion::ComputeTime ( double  sTransform,
double  D,
double  r 
)
inlinestatic

Definition at line 87 of file G4DNASmoluchowskiDiffusion.hh.

88  {
89  return std::pow(r / sTransform, 2.) / (4. * D);
90  }
double D(double temp)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 110 of file G4DNASmoluchowskiDiffusion.hh.

113  {
114  double sTransform = GetInverseProbability(proba);
115  return ComputeTime(sTransform, D, distance);
116  }
static double ComputeTime(double sTransform, double D, double r)
double D(double temp)

Here is the call graph for this function:

static double G4DNASmoluchowskiDiffusion::GetCumulativeProbability ( double  sTransform)
inlinestatic

Definition at line 304 of file G4DNASmoluchowskiDiffusion.hh.

305  {
306  static double constant = 2./std::sqrt(3.141592653589793);
307  return erfc(sTransform) + constant*sTransform*G4Exp(-sTransform*sTransform);
308  }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 129 of file G4DNASmoluchowskiDiffusion.hh.

130  {
131  static double my_pi = 3.141592653589793;
132  static double constant = 4.*my_pi/std::pow(4.*my_pi, 1.5);
133  return r*r/std::pow(D * _time,1.5)*G4Exp(-r*r/(4. * D * _time))*constant;
134  }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double D(double temp)

Here is the call graph for this function:

static double G4DNASmoluchowskiDiffusion::GetDifferential ( double  sTransform)
inlinestatic

Definition at line 123 of file G4DNASmoluchowskiDiffusion.hh.

124  {
125  static double constant = -4./std::sqrt(3.141592653589793);
126  return sTransform*sTransform*G4Exp(-sTransform*sTransform)*constant; // -4*sTransform*sTransform*exp(-sTransform*sTransform)/sqrt(3.141592653589793)
127  }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183

Here is the call graph for this function:

Here is the caller graph for this function:

double G4DNASmoluchowskiDiffusion::GetInverseProbability ( double  proba)
inline

Definition at line 310 of file G4DNASmoluchowskiDiffusion.hh.

311  {
312  size_t index_low = (size_t) trunc(proba/fEpsilon);
313 
314  if(index_low == 0) // assymptote en 0
315  {
316  index_low = 1;
317  size_t index_up = 2;
318  double low_y = fInverse[index_low];
319  double up_y = fInverse[index_up];
320  double low_x = index_low*fEpsilon;
321  double up_x = proba+fEpsilon;
322  double tangente = (low_y-up_y)/(low_x - up_x); // ou utiliser GetDifferential(proba) ?
323  // double tangente = GetDifferential(proba);
324  return low_y + tangente*(proba-low_x);
325  }
326 
327  size_t index_up = index_low+1;
328  if(index_low > fInverse.size()) return fInverse.back();
329  double low_y = fInverse[index_low];
330  double up_y = fInverse[index_up];
331 
332  double low_x = index_low*fEpsilon;
333  double up_x = low_x+fEpsilon;
334 
335  if(up_x > 1) // P(1) = 0
336  {
337  up_x = 1;
338  up_y = 0; // more general : fInverse.back()
339  }
340 
341  double tangente = (low_y-up_y)/(low_x - up_x);
342 
343  return low_y + tangente*(proba-low_x);
344  }

Here is the caller graph for this function:

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

Definition at line 94 of file G4DNASmoluchowskiDiffusion.hh.

95  {
96  double proba = G4UniformRand();
97 // G4cout << "proba = " << proba << G4endl;
98  double sTransform = GetInverseProbability(proba);
99 // G4cout << "sTransform = " << sTransform << G4endl;
100  return ComputeDistance(sTransform, D, _time);
101  }
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:

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

Definition at line 103 of file G4DNASmoluchowskiDiffusion.hh.

104  {
105  double proba = G4UniformRand();
106  double sTransform = GetInverseProbability(proba);
107  return ComputeTime(sTransform, D, distance);
108  }
#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:

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

Definition at line 357 of file G4DNASmoluchowskiDiffusion.hh.

358  {
359  // x > x'
360  // P'(x) = p(x') = lim(x->x') (P(x') - P(x))/(x'-x)
361  // x'-x = (P(x') - P(x))/p(x')
362  // x = x' - (P(x') - P(x))/p(x')
363 
364  // fInverse initialized in the constructor
365 
366  assert(fNbins !=0);
367  PrepareReverseTable(0,xmax);
368  }
void PrepareReverseTable(double xmin, double xmax)

Here is the call graph for this function:

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

Definition at line 351 of file G4DNASmoluchowskiDiffusion.hh.

352  {
353  return GetDifferential(x[0]);
354  }
static double GetDifferential(double sTransform)

Here is the call graph for this function:

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

Definition at line 346 of file G4DNASmoluchowskiDiffusion.hh.

347  {
348  return GetInverseProbability(x[0]);
349  }

Here is the call graph for this function:

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

Definition at line 263 of file G4DNASmoluchowskiDiffusion.hh.

264  {
265  double x = xmax;
266  int index = 0;
267  double nextProba = fEpsilon;
268  double proposedX;
269 
270  BoundingBox boundingBox(xmin, xmax, fEpsilon*1e-5);
271 
272  while(index <= fNbins)
273  // in case GetCumulativeProbability is exact (digitally speaking), replace with:
274  // while(index <= fNbins+1)
275  {
276  nextProba = fEpsilon*index;
277 
278  double newProba = GetCumulativeProbability(x);
279 
280  if(boundingBox.Propose(x, newProba, nextProba, proposedX))
281  {
282  fInverse[index] = x;
283  index++;
284  }
285  else
286  {
287  if(x == proposedX)
288  {
289  G4cout << "BREAK : x= " << x << G4endl;
290  G4cout << "index= " << index << G4endl;
291  G4cout << "nextProba= " << nextProba << G4endl;
292  G4cout << "newProba= " << newProba << G4endl;
293  abort();
294  }
295  x = proposedX;
296  }
297  }
298 
299  fInverse[fNbins+1] = 0; // P(1) = 0, because we want it exact !
300  // Tips to improve the exactness: get an better value of pi, get better approximation of erf and exp, use long double instead of double
301 // boundingBox.Print();
302  }
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

double G4DNASmoluchowskiDiffusion::fEpsilon

Definition at line 372 of file G4DNASmoluchowskiDiffusion.hh.

std::vector<double> G4DNASmoluchowskiDiffusion::fInverse

Definition at line 370 of file G4DNASmoluchowskiDiffusion.hh.

int G4DNASmoluchowskiDiffusion::fNbins

Definition at line 371 of file G4DNASmoluchowskiDiffusion.hh.


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