Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingFastSimInterpolation.cc
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/// \file G4ChannelingFastSimInterpolation.cc
28/// \brief Implementation of the G4ChannelingFastSimInterpolation class
29
31#include "G4SystemOfUnits.hh"
32
34 G4double dy0,
35 G4int nPointsx0,
36 G4int nPointsy0,
37 G4int iModel0)
38{
39 fDx = dx0;
40 fDy = dy0;
41 nPointsx = nPointsx0;
42 nPointsy = nPointsy0;
43 fStepi = fDx/nPointsx;
44 fStepi2= fStepi*fStepi;
45 iModel = iModel0;
46
47 //resize vectors for interpolation coefficients
48 if(iModel==1)
49 {
50 fAI.resize(nPointsx+1);
51 fBI.resize(nPointsx+1);
52 fCI.resize(nPointsx+1);
53 fDI.resize(nPointsx+1);
54 }
55 else if(iModel==2)
56 {
57 fStepj = fDy/nPointsy;
58 fAI3D.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
59 fBI3D.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
60 fCI3D.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
61 fAI3D3.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
62 fBI3D3.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
63 fCI3D3.resize(nPointsx+1, std::vector<G4double>(nPointsy+1));
64
65 for(G4int i=0; i<=nPointsx; i++)
66 {
67 fCI3D[i][nPointsy]=0.;
68 }
69 }
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75 G4double SplineF=0.;
76 if(iModel==1)
77 {
78 SplineF = Spline1D(xx);
79 }
80 else if(iModel==2)
81 {
82 SplineF = Spline2D(xx, yy);
83 }
84 return SplineF;
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90 G4double BI0,
91 G4double CI0,
92 G4double DI0,
93 G4int i){
94 fAI[i] = AI0;
95 fBI[i] = BI0/cm;
96 fCI[i] = CI0/cm2;
97 fDI[i] = DI0/cm3;
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
103 G4double BI3D0,
104 G4double CI3D0,
105 G4int i,
106 G4int j,
107 G4int k){
108 if (k==0)
109 {
110 fAI3D[i][j] = AI3D0/fStepj/fStepi/6.;
111 fBI3D[i][j] = BI3D0/fStepj/fStepi/6.;
112 fCI3D[i][j] = CI3D0/fStepj/fStepi/6./cm2;
113 }
114 else if (k==1)
115 {
116 fAI3D3[i][j] = AI3D0/fStepj/fStepi/6./cm2;
117 fBI3D3[i][j] = BI3D0/fStepj/fStepi/6./cm2;
118 fCI3D3[i][j] = CI3D0/fStepj/fStepi/6./cm2/cm2;
119 }
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
124G4double G4ChannelingFastSimInterpolation::Spline1D(G4double xx) //cubic spline of
125 //1-variable function
126{
127 G4double x1 = xx;
128
129 //if a particle escapes the interpolation area
130 if (x1<0.)
131 {
132 x1 += fDx;
133 }
134 else if (x1>=fDx)
135 {
136 x1 -= fDx;
137 }
138
139// calculation of interpolation function
140 G4int mmx = std::floor(x1/fStepi);
141 x1 -= (mmx+1)*fStepi;
142 G4double Spline3 = fAI[mmx]+x1*(fBI[mmx]+x1*(fCI[mmx]+fDI[mmx]*x1));
143 return Spline3;
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
148 G4double G4ChannelingFastSimInterpolation::Spline2D(G4double xx, G4double yy)//cubic
149 //spline of 2-variable function
150 {
151 G4double x1 = xx;
152 G4double y1 = yy;
153
154 //if a particle escapes the interpolation area
155 if (x1<0.)
156 {
157 x1 += fDx;
158 }
159 else if (x1>=fDx)
160 {
161 x1 -= fDx;
162 }
163 if (y1<0.)
164 {
165 y1 += fDy;
166 }
167 else if (y1>=fDy)
168 {
169 y1 -= fDy;
170 }
171
172 // calculation of interpolation function
173 G4int mmx = std::floor(x1/fStepi);
174 G4int mmy = std::floor(y1/fStepj);
175
176 G4double tt1 = x1-mmx*fStepi;
177 G4double tt2 = fStepi-tt1;
178 G4double tt13 = tt1*tt1*tt1;
179 G4double tt23 = tt2*tt2*tt2;
180
181 G4double tt1y = y1-mmy*fStepj;
182 G4double tt2y = fStepj-tt1y;
183 G4double tt1y3 = tt1y*tt1y*tt1y;
184 G4double tt2y3 = tt2y*tt2y*tt2y;
185
186 G4double spl3dxx1 = fCI3D3[mmx ][mmy]*tt2y3 + fCI3D3[mmx ][mmy+1]*tt1y3 +
187 fAI3D3[mmx ][mmy]*tt2y + fBI3D3[mmx ][mmy]*tt1y;
188 G4double spl3dxx2 = fCI3D3[mmx+1][mmy]*tt2y3 + fCI3D3[mmx+1][mmy+1]*tt1y3 +
189 fAI3D3[mmx+1][mmy]*tt2y + fBI3D3[mmx+1][mmy]*tt1y;
190 G4double spl3d1 = fCI3D[ mmx ][mmy]*tt2y3 + fCI3D[mmx ][mmy+1]*tt1y3 +
191 fAI3D[ mmx ][mmy]*tt2y + fBI3D[mmx ][mmy]*tt1y;
192 G4double spl3d2 = fCI3D[ mmx+1][mmy]*tt2y3 + fCI3D[mmx+1][mmy+1]*tt1y3 +
193 fAI3D[ mmx+1][mmy]*tt2y + fBI3D[mmx+1][mmy]*tt1y;
194
195 G4double spl3d = spl3dxx1*tt23 + spl3dxx2*tt13 +
196 (spl3d1*6.-spl3dxx1*fStepi2)*tt2 + (spl3d2*6.-spl3dxx2*fStepi2)*tt1;
197
198 return spl3d;
199 }
200
201
Definition of the G4ChannelingFastSimInterpolation class The class includes spline interpolation coef...
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void SetCoefficients1D(G4double AI0, G4double BI0, G4double CI0, G4double DI0, G4int i)
Set spline coefficients.
void SetCoefficients2D(G4double AI3D0, G4double BI3D0, G4double CI3D0, G4int i, G4int j, G4int k)
G4ChannelingFastSimInterpolation(G4double dx0, G4double dy0, G4int nPointsx0, G4int nPointsy0, G4int iModel0)
G4double GetIF(G4double xx, G4double yy)
Get Spline Function.