71double time_test = 1e-6 ;
73double test_distance = 1e-9;
75double Plot(
double* x,
double* )
81static double InvErfc(
double x)
83 return TMath::ErfcInverse(x);
86Axis_t* BinLogX(Int_t bins, Axis_t from, Axis_t to)
88 Axis_t width = (to - from) / bins;
89 Axis_t *new_bins =
new Axis_t[bins + 1];
91 for (
int i = 0; i <= bins; i++) {
92 new_bins[i] = TMath::Power(10, from + i * width);
98int main(
int argc,
char **argv)
102 TRint* root =
new TRint(
"G4DNASmoluchowskiDiffusion",&argc, argv);
103 double interval = 1e-5;
112 std::cout << diff->
fInverse.size() << std::endl;
114 TCanvas* canvas =
new TCanvas();
124 TH1D* h1 =
new TH1D(
"h1",
"h1", 100, 0., 1e-6);
125 double distance = -1;
129 for(
size_t i = 0 ; i <
N ; ++i)
139 int integral_h1 = h1->Integral();
140 h1->Scale(1./integral_h1);
141 scalf=h1->GetBinWidth ( 1 ) ;
143 h1->GetXaxis()->SetTitle(
"distance");
146 TH1D* h2 =
new TH1D(
"h2",
"h2", 100, 0., 1e-6);
147 TH1D* h_irt_distance =
new TH1D(
"h2",
"h2", 100, 0., 1e-6);
149 for(
size_t i = 0 ; i <
N ; ++i)
151 double x = std::sqrt(2*
D*time_test)*root_random.Gaus();
152 double y = std::sqrt(2*
D*time_test)*root_random.Gaus();
153 double z = std::sqrt(2*
D*time_test)*root_random.Gaus();
155 distance = std::sqrt(x*x+y*y+z*z);
159 double proba = root_random.Rndm();
160 double irt_distance = InvErfc(proba)*2*std::sqrt(
D*time_test);
161 h_irt_distance->Fill(irt_distance);
165 int integral_h2 = h2->Integral();
166 h2->Scale(1./integral_h2);
167 scalf=h2->GetBinWidth ( 1 ) ;
172 int integral_h_irt_distance = h_irt_distance->Integral();
173 h_irt_distance->Scale(1./integral_h_irt_distance);
174 scalf = h_irt_distance->GetBinWidth ( 1 ) ;
175 h_irt_distance->Scale(1./scalf);
176 h_irt_distance->GetXaxis()->SetTitle(
"distance");
180 TF1 * f2 =
new TF1(
"f2",&Plot,0,1e-6,0,
"Plot");
186 h_irt_distance->Draw(
"SAME");
187 double integral = f2->Integral(0., 1e-6);
188 std::cout <<
"integral = " << integral << std::endl;
189 std::cout <<
"integral h1 = " << h1->Integral() << std::endl;
192 std::vector<double> rdm(3);
194 Axis_t* bins = BinLogX(nbins, -12, -1);
196 TH1D* h3 =
new TH1D(
"h3",
"h3", 100, bins);
197 TH1D* h4 =
new TH1D(
"h4",
"h4", 100, bins);
198 TH1D* h_irt =
new TH1D(
"h_irt",
"h_irt", 100, bins);
200 for(
size_t i = 0 ; i <
N ; ++i)
202 for(
size_t j = 0 ; j < 3 ; ++j)
203 rdm[j] = root_random.Gaus();
205 double denum = 1./(rdm[0]*rdm[0] + rdm[1]*rdm[1] + rdm[2]*rdm[2]);
207 double t = ((test_distance*test_distance)*denum)*1./(2*
D);
214 double proba = root_random.Rndm();
215 double t_irt = 1./(4*
D)*std::pow((test_distance)/InvErfc(proba),2);
220 TCanvas* c1 =
new TCanvas();
222 int integral_h3 = h3->Integral();
223 h3->Scale(1./integral_h3);
224 scalf=h3->GetBinWidth ( 1 ) ;
227 h3->GetXaxis()->SetTitle(
"time");;
234 int integral_h4 = h4->Integral();
235 h4->Scale(1./integral_h4);
236 scalf=h4->GetBinWidth ( 1 ) ;
246 int integral_h_irt = h_irt->Integral();
247 h_irt->Scale(1./integral_h_irt);
248 scalf=h_irt->GetBinWidth ( 1 ) ;
249 h_irt->Scale(1./scalf);
250 h_irt->SetLineColor(4);
G4double epsilon(G4double density, G4double temperature)
G4double D(G4double temp)
std::vector< double > fInverse
static double GetDensityProbability(double r, double _time, double D)
double GetRandomTime(double distance, double D)
virtual ~G4DNASmoluchowskiDiffusion()
void InitialiseInverseProbability(double xmax=3e28)
G4DNASmoluchowskiDiffusion(double epsilon=1e-5)
double GetRandomDistance(double _time, double D)