Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPLabAngularEnergy Class Reference

#include <G4ParticleHPLabAngularEnergy.hh>

+ Inheritance diagram for G4ParticleHPLabAngularEnergy:

Public Member Functions

 G4ParticleHPLabAngularEnergy ()
 
 ~G4ParticleHPLabAngularEnergy () override
 
void Init (std::istream &aDataFile) override
 
G4ReactionProductSample (G4double anEnergy, G4double massCode, G4double mass) override
 
G4double MeanEnergyOfThisInteraction () override
 
- Public Member Functions inherited from G4VParticleHPEnergyAngular
 G4VParticleHPEnergyAngular ()
 
virtual ~G4VParticleHPEnergyAngular ()=default
 
void SetProjectileRP (G4ReactionProduct *aIncidentParticleRP)
 
void SetTarget (G4ReactionProduct *aTarget)
 
G4ReactionProductGetTarget () const
 
G4ReactionProductGetProjectileRP () const
 
G4ReactionProductGetCMS ()
 
void SetQValue (G4double aValue)
 
virtual void ClearHistories ()
 
 G4VParticleHPEnergyAngular (G4VParticleHPEnergyAngular &)=delete
 
G4VParticleHPEnergyAngularoperator= (const G4VParticleHPEnergyAngular &right)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleHPEnergyAngular
G4double GetQValue () const
 

Detailed Description

Definition at line 43 of file G4ParticleHPLabAngularEnergy.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPLabAngularEnergy()

G4ParticleHPLabAngularEnergy::G4ParticleHPLabAngularEnergy ( )
inline

Definition at line 46 of file G4ParticleHPLabAngularEnergy.hh.

47 {
48 theEnergies = nullptr;
49 theData = nullptr;
50 nCosTh = nullptr;
51 theSecondManager = nullptr;
52 nEnergies = -1;
53 currentMeanEnergy = -1.0;
54 }

◆ ~G4ParticleHPLabAngularEnergy()

G4ParticleHPLabAngularEnergy::~G4ParticleHPLabAngularEnergy ( )
inlineoverride

Definition at line 55 of file G4ParticleHPLabAngularEnergy.hh.

56 {
57 delete[] theEnergies;
58 delete[] nCosTh;
59 if (theData != nullptr) {
60 for (G4int i = 0; i < nEnergies; i++)
61 delete[] theData[i];
62 delete[] theData;
63 }
64 delete[] theSecondManager;
65 }
int G4int
Definition G4Types.hh:85

Member Function Documentation

◆ Init()

void G4ParticleHPLabAngularEnergy::Init ( std::istream & aDataFile)
overridevirtual

Implements G4VParticleHPEnergyAngular.

Definition at line 51 of file G4ParticleHPLabAngularEnergy.cc.

52{
53 aDataFile >> nEnergies;
54 theManager.Init(aDataFile);
55 const std::size_t esize = nEnergies > 0 ? nEnergies : 1;
56 theEnergies = new G4double[esize];
57 nCosTh = new G4int[esize];
58 theData = new G4ParticleHPVector*[esize];
59 theSecondManager = new G4InterpolationManager[esize];
60 for (G4int i = 0; i < nEnergies; ++i) {
61 aDataFile >> theEnergies[i];
62 theEnergies[i] *= eV;
63 aDataFile >> nCosTh[i];
64 theSecondManager[i].Init(aDataFile);
65 const std::size_t dsize = nCosTh[i] > 0 ? nCosTh[i] : 1;
66 theData[i] = new G4ParticleHPVector[dsize];
67 G4double label;
68 for (std::size_t ii = 0; ii < dsize; ++ii) {
69 aDataFile >> label;
70 theData[i][ii].SetLabel(label);
71 theData[i][ii].Init(aDataFile, eV);
72 }
73 }
74}
double G4double
Definition G4Types.hh:83
void Init(G4int aScheme, G4int aRange)
void SetLabel(G4double aLabel)
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)

◆ MeanEnergyOfThisInteraction()

G4double G4ParticleHPLabAngularEnergy::MeanEnergyOfThisInteraction ( )
inlineoverridevirtual

Implements G4VParticleHPEnergyAngular.

Definition at line 70 of file G4ParticleHPLabAngularEnergy.hh.

70{ return currentMeanEnergy; }

◆ Sample()

G4ReactionProduct * G4ParticleHPLabAngularEnergy::Sample ( G4double anEnergy,
G4double massCode,
G4double mass )
overridevirtual

Implements G4VParticleHPEnergyAngular.

Definition at line 76 of file G4ParticleHPLabAngularEnergy.cc.

78{
79 auto result = new G4ReactionProduct;
80 auto Z = static_cast<G4int>(massCode / 1000);
81 auto A = static_cast<G4int>(massCode - 1000 * Z);
82
83 if (massCode == 0) {
84 result->SetDefinition(G4Gamma::Gamma());
85 }
86 else if (A == 0) {
87 result->SetDefinition(G4Electron::Electron());
88 if (Z == 1) result->SetDefinition(G4Positron::Positron());
89 }
90 else if (A == 1) {
91 result->SetDefinition(G4Neutron::Neutron());
92 if (Z == 1) result->SetDefinition(G4Proton::Proton());
93 }
94 else if (A == 2) {
95 result->SetDefinition(G4Deuteron::Deuteron());
96 }
97 else if (A == 3) {
98 result->SetDefinition(G4Triton::Triton());
99 if (Z == 2) result->SetDefinition(G4He3::He3());
100 }
101 else if (A == 4) {
102 result->SetDefinition(G4Alpha::Alpha());
103 if (Z != 2) throw G4HadronicException(__FILE__, __LINE__, "Unknown ion case 1");
104 }
105 else {
106 throw G4HadronicException(__FILE__, __LINE__,
107 "G4ParticleHPLabAngularEnergy: Unknown ion case 2");
108 }
109
110 // get theta, E
111 G4double cosTh, secEnergy;
112 G4int i, it(0);
113 // find the energy bin
114 for (i = 0; i < nEnergies; i++) {
115 it = i;
116 if (anEnergy < theEnergies[i]) break;
117 }
118
119 if (it == 0) // it marks the energy bin --> we do not extrapolate to low energies, we extrapolate
120 // to high energies (??)
121 {
122 // Do not sample between incident neutron energies:
123 //---------------------------------------------------------
124 // CosTheta:
125 G4ParticleHPVector theThVec;
126 theThVec.SetInterpolationManager(theSecondManager[it]);
127 for (i = 0; i < nCosTh[it]; i++) {
128 theThVec.SetX(i, theData[it][i].GetLabel());
129 theThVec.SetY(i, theData[it][i].GetIntegral());
130 }
131 cosTh = theThVec.Sample();
132 //---------------------------------------------------------
133
134 //---------------------------------------------------------
135 // Now the secondary energy:
136 G4double x, x1, x2, y1, y2, y, E;
137 G4int i1(0);
138 for (i = 1; i < nCosTh[it]; i++) {
139 i1 = i;
140 if (cosTh < theData[it][i].GetLabel()) break;
141 }
142 // now get the prob at this energy for the right theta value
143 x = cosTh;
144 x1 = theData[it][i1 - 1].GetLabel();
145 x2 = theData[it][i1].GetLabel();
146 G4ParticleHPVector theBuff1a;
147 theBuff1a.SetInterpolationManager(theData[it][i1 - 1].GetInterpolationManager());
148 for (i = 0; i < theData[it][i1 - 1].GetVectorLength(); i++) {
149 E = theData[it][i1 - 1].GetX(i);
150 y1 = theData[it][i1 - 1].GetY(i);
151 y2 = theData[it][i1].GetY(E);
152 y = theInt.Interpolate(theSecondManager[it].GetScheme(i1), x, x1, x2, y1, y2);
153 theBuff1a.SetData(i, E, y);
154 }
155 G4ParticleHPVector theBuff2a;
156 theBuff2a.SetInterpolationManager(theData[it][i1].GetInterpolationManager());
157 for (i = 0; i < theData[it][i1].GetVectorLength(); i++) {
158 E = theData[it][i1].GetX(i);
159 y1 = theData[it][i1 - 1].GetY(E);
160 y2 = theData[it][i1].GetY(i);
161 y = theInt.Interpolate(theSecondManager[it].GetScheme(i1), x, x1, x2, y1, y2);
162 theBuff2a.SetData(i, E, y); // wrong E, right theta.
163 }
164 G4ParticleHPVector theStore;
165 theStore.Merge(&theBuff1a, &theBuff2a);
166 secEnergy = theStore.Sample();
167 currentMeanEnergy = theStore.GetMeanX();
168 //---------------------------------------------------------
169 }
170 else // this is the small big else.
171 {
172 G4double x, x1, x2, y1, y2, y, tmp, E;
173 // integrate the prob for each costh, and select theta.
175 for (i = 0; i < nCosTh[it - 1]; i++) {
176 run1.SetX(i, theData[it - 1][i].GetLabel());
177 run1.SetY(i, theData[it - 1][i].GetIntegral());
178 }
180 for (i = 0; i < nCosTh[it]; i++) {
181 run2.SetX(i, theData[it][i].GetLabel());
182 run2.SetY(i, theData[it][i].GetIntegral());
183 }
184 // get the distributions for the correct neutron energy
185 x = anEnergy;
186 x1 = theEnergies[it - 1];
187 x2 = theEnergies[it];
188 G4ParticleHPVector thBuff1; // to be interpolated as run1.
189 thBuff1.SetInterpolationManager(theSecondManager[it - 1]);
190 for (i = 0; i < run1.GetVectorLength(); i++) {
191 tmp = run1.GetX(i); // theta
192 y1 = run1.GetY(i); // integral
193 y2 = run2.GetY(tmp);
194 y = theInt.Interpolate(theManager.GetScheme(it), x, x1, x2, y1, y2);
195 thBuff1.SetData(i, tmp, y);
196 }
197 G4ParticleHPVector thBuff2;
198 thBuff2.SetInterpolationManager(theSecondManager[it]);
199 for (i = 0; i < run2.GetVectorLength(); i++) {
200 tmp = run2.GetX(i); // theta
201 y1 = run1.GetY(tmp); // integral
202 y2 = run2.GetY(i);
203 y = theInt.Interpolate(theManager.GetScheme(it), x, x1, x2, y1, y2);
204 thBuff2.SetData(i, tmp, y);
205 }
206 G4ParticleHPVector theThVec;
207 theThVec.Merge(&thBuff1, &thBuff2); // takes care of interpolation
208 cosTh = theThVec.Sample();
209 /*
210 if(massCode>0.99 && massCode<1.01){theThVec.Dump();}
211 G4double random = (theThVec.GetY(theThVec.GetVectorLength()-1)
212 -theThVec.GetY(0)) *G4UniformRand();
213 G4cout<<" -- "<<theThVec.GetY(theThVec.GetVectorLength()-1)<<" "<<theThVec.GetY(0)<<" ---->
214 "<<random<<G4endl; G4int ith(0); for(i=1;i<theThVec.GetVectorLength(); i++)
215 {
216 ith = i;
217 if(random<theThVec.GetY(i)-theThVec.GetY(0)) break;
218 }
219 {
220 // calculate theta
221 G4double xx, xx1, xx2, yy1, yy2;
222 xx = random;
223 xx1 = theThVec.GetY(ith-1)-theThVec.GetY(0); // integrals
224 xx2 = theThVec.GetY(ith)-theThVec.GetY(0);
225 yy1 = theThVec.GetX(ith-1); // std::cos(theta)
226 yy2 = theThVec.GetX(ith);
227 cosTh = theInt.Interpolate(theSecondManager[it].GetScheme(ith),
228 xx, xx1,xx2,yy1,yy2);
229 }
230 */
231 G4int i1(0), i2(0);
232 // get the indixes of the vectors close to theta for low energy
233 // first it-1 !!!! i.e. low in energy
234 for (i = 1; i < nCosTh[it - 1]; i++) {
235 i1 = i;
236 if (cosTh < theData[it - 1][i].GetLabel()) break;
237 }
238 // now get the prob at this energy for the right theta value
239 x = cosTh;
240 x1 = theData[it - 1][i1 - 1].GetLabel();
241 x2 = theData[it - 1][i1].GetLabel();
242 G4ParticleHPVector theBuff1a;
243 theBuff1a.SetInterpolationManager(theData[it - 1][i1 - 1].GetInterpolationManager());
244 for (i = 0; i < theData[it - 1][i1 - 1].GetVectorLength(); i++) {
245 E = theData[it - 1][i1 - 1].GetX(i);
246 y1 = theData[it - 1][i1 - 1].GetY(i);
247 y2 = theData[it - 1][i1].GetY(E);
248 y = theInt.Interpolate(theSecondManager[it - 1].GetScheme(i1), x, x1, x2, y1, y2);
249 theBuff1a.SetData(i, E, y); // wrong E, right theta.
250 }
251 G4ParticleHPVector theBuff2a;
252 theBuff2a.SetInterpolationManager(theData[it - 1][i1].GetInterpolationManager());
253 for (i = 0; i < theData[it - 1][i1].GetVectorLength(); i++) {
254 E = theData[it - 1][i1].GetX(i);
255 y1 = theData[it - 1][i1 - 1].GetY(E);
256 y2 = theData[it - 1][i1].GetY(i);
257 y = theInt.Interpolate(theSecondManager[it - 1].GetScheme(i1), x, x1, x2, y1, y2);
258 theBuff2a.SetData(i, E, y); // wrong E, right theta.
259 }
260 G4ParticleHPVector theStore1;
261 theStore1.Merge(&theBuff1a, &theBuff2a); // wrong E, right theta, complete binning
262
263 // get the indixes of the vectors close to theta for high energy
264 // then it !!!! i.e. high in energy
265 for (i = 1; i < nCosTh[it]; i++) {
266 i2 = i;
267 if (cosTh < theData[it][i2].GetLabel()) break;
268 } // sonderfaelle mit i1 oder i2 head on fehlen. @@@@@
269 x1 = theData[it][i2 - 1].GetLabel();
270 x2 = theData[it][i2].GetLabel();
271 G4ParticleHPVector theBuff1b;
272 theBuff1b.SetInterpolationManager(theData[it][i2 - 1].GetInterpolationManager());
273 for (i = 0; i < theData[it][i2 - 1].GetVectorLength(); i++) {
274 E = theData[it][i2 - 1].GetX(i);
275 y1 = theData[it][i2 - 1].GetY(i);
276 y2 = theData[it][i2].GetY(E);
277 y = theInt.Interpolate(theSecondManager[it].GetScheme(i2), x, x1, x2, y1, y2);
278 theBuff1b.SetData(i, E, y); // wrong E, right theta.
279 }
280 G4ParticleHPVector theBuff2b;
281 theBuff2b.SetInterpolationManager(theData[it][i2].GetInterpolationManager());
282 // 080808 i1 -> i2
283 // for(i=0;i<theData[it][i1].GetVectorLength(); i++)
284 for (i = 0; i < theData[it][i2].GetVectorLength(); i++) {
285 // E = theData[it][i1].GetX(i);
286 // y1 = theData[it][i1-1].GetY(E);
287 // y2 = theData[it][i1].GetY(i);
288 E = theData[it][i2].GetX(i);
289 y1 = theData[it][i2 - 1].GetY(E);
290 y2 = theData[it][i2].GetY(i);
291 y = theInt.Interpolate(theSecondManager[it].GetScheme(i2), x, x1, x2, y1, y2);
292 theBuff2b.SetData(i, E, y); // wrong E, right theta.
293 }
294 G4ParticleHPVector theStore2;
295 theStore2.Merge(&theBuff1b, &theBuff2b); // wrong E, right theta, complete binning
296 // now get to the right energy.
297
298 x = anEnergy;
299 x1 = theEnergies[it - 1];
300 x2 = theEnergies[it];
301 G4ParticleHPVector theOne1;
303 for (i = 0; i < theStore1.GetVectorLength(); i++) {
304 E = theStore1.GetX(i);
305 y1 = theStore1.GetY(i);
306 y2 = theStore2.GetY(E);
307 y = theInt.Interpolate(theManager.GetScheme(it), x, x1, x2, y1, y2);
308 theOne1.SetData(i, E, y); // both correct
309 }
310 G4ParticleHPVector theOne2;
312 for (i = 0; i < theStore2.GetVectorLength(); i++) {
313 E = theStore2.GetX(i);
314 y1 = theStore1.GetY(E);
315 y2 = theStore2.GetY(i);
316 y = theInt.Interpolate(theManager.GetScheme(it), x, x1, x2, y1, y2);
317 theOne2.SetData(i, E, y); // both correct
318 }
319 G4ParticleHPVector theOne;
320 theOne.Merge(&theOne1, &theOne2); // both correct, complete binning
321
322 secEnergy = theOne.Sample();
323 currentMeanEnergy = theOne.GetMeanX();
324 }
325
326 // now do random direction in phi, and fill the result.
327
328 result->SetKineticEnergy(secEnergy);
329
330 G4double phi = twopi * G4UniformRand();
331 G4double theta = std::acos(cosTh);
332 G4double sinth = std::sin(theta);
333 G4double mtot = result->GetTotalMomentum();
334 G4ThreeVector tempVector(mtot * sinth * std::cos(phi), mtot * sinth * std::sin(phi),
335 mtot * std::cos(theta));
336 result->SetMomentum(tempVector);
337
338 return result;
339}
const G4double A[17]
#define G4UniformRand()
Definition Randomize.hh:52
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4He3 * He3()
Definition G4He3.cc:90
G4InterpolationScheme GetScheme(G4int index) const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4double Interpolate(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void SetY(G4int i, G4double x)
void SetX(G4int i, G4double e)
void SetData(G4int i, G4double x, G4double y)
const G4InterpolationManager & GetInterpolationManager() const
void SetInterpolationManager(const G4InterpolationManager &aManager)
G4double GetY(G4double x)
G4double GetX(G4int i) const
void Merge(G4ParticleHPVector *active, G4ParticleHPVector *passive)
G4int GetVectorLength() const
static G4Positron * Positron()
Definition G4Positron.cc:90
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4Triton * Triton()
Definition G4Triton.cc:90

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