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 | // $Id: G4MuBremsstrahlungModel.cc,v 1.35 2009/04/12 17:48:45 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4MuBremsstrahlungModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
37 | // |
---|
38 | // Creation date: 24.06.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko) |
---|
43 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) |
---|
44 | // 24-01-03 Fix for compounds (V.Ivanchenko) |
---|
45 | // 27-01-03 Make models region aware (V.Ivanchenko) |
---|
46 | // 13-02-03 Add name (V.Ivanchenko) |
---|
47 | // 10-02-04 Add lowestKinEnergy (V.Ivanchenko) |
---|
48 | // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko) |
---|
49 | // 03-08-05 Angular correlations according to PRM (V.Ivanchenko) |
---|
50 | // 13-02-06 add ComputeCrossSectionPerAtom (mma) |
---|
51 | // 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI) |
---|
52 | // 07-11-07 Improve sampling of final state (A.Bogdanov) |
---|
53 | // 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko) |
---|
54 | // |
---|
55 | |
---|
56 | // |
---|
57 | // Class Description: |
---|
58 | // |
---|
59 | // |
---|
60 | // ------------------------------------------------------------------- |
---|
61 | // |
---|
62 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
64 | |
---|
65 | #include "G4MuBremsstrahlungModel.hh" |
---|
66 | #include "G4Gamma.hh" |
---|
67 | #include "G4MuonMinus.hh" |
---|
68 | #include "G4MuonPlus.hh" |
---|
69 | #include "Randomize.hh" |
---|
70 | #include "G4Material.hh" |
---|
71 | #include "G4Element.hh" |
---|
72 | #include "G4ElementVector.hh" |
---|
73 | #include "G4ProductionCutsTable.hh" |
---|
74 | #include "G4ParticleChangeForLoss.hh" |
---|
75 | |
---|
76 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
77 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
78 | |
---|
79 | using namespace std; |
---|
80 | |
---|
81 | G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(const G4ParticleDefinition* p, |
---|
82 | const G4String& nam) |
---|
83 | : G4VEmModel(nam), |
---|
84 | particle(0), |
---|
85 | sqrte(sqrt(exp(1.))), |
---|
86 | bh(202.4), |
---|
87 | bh1(446.), |
---|
88 | btf(183.), |
---|
89 | btf1(1429.), |
---|
90 | fParticleChange(0), |
---|
91 | lowestKinEnergy(1.0*GeV), |
---|
92 | minThreshold(1.0*keV) |
---|
93 | { |
---|
94 | theGamma = G4Gamma::Gamma(); |
---|
95 | nist = G4NistManager::Instance(); |
---|
96 | if(p) SetParticle(p); |
---|
97 | } |
---|
98 | |
---|
99 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
100 | |
---|
101 | G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel() |
---|
102 | { |
---|
103 | size_t n = partialSumSigma.size(); |
---|
104 | if(n > 0) { |
---|
105 | for(size_t i=0; i<n; i++) { |
---|
106 | delete partialSumSigma[i]; |
---|
107 | } |
---|
108 | } |
---|
109 | } |
---|
110 | |
---|
111 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
112 | |
---|
113 | G4double G4MuBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*, |
---|
114 | const G4MaterialCutsCouple*) |
---|
115 | { |
---|
116 | return minThreshold; |
---|
117 | } |
---|
118 | |
---|
119 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
120 | |
---|
121 | void G4MuBremsstrahlungModel::Initialise(const G4ParticleDefinition* p, |
---|
122 | const G4DataVector& cuts) |
---|
123 | { |
---|
124 | if(p) SetParticle(p); |
---|
125 | |
---|
126 | highKinEnergy = HighEnergyLimit(); |
---|
127 | |
---|
128 | // partial cross section is computed for fixed energy |
---|
129 | G4double fixedEnergy = 0.5*highKinEnergy; |
---|
130 | |
---|
131 | const G4ProductionCutsTable* theCoupleTable= |
---|
132 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
133 | if(theCoupleTable) { |
---|
134 | G4int numOfCouples = theCoupleTable->GetTableSize(); |
---|
135 | |
---|
136 | // clear old data |
---|
137 | G4int nn = partialSumSigma.size(); |
---|
138 | G4int nc = cuts.size(); |
---|
139 | if(nn > 0) { |
---|
140 | for (G4int ii=0; ii<nn; ii++){ |
---|
141 | G4DataVector* a = partialSumSigma[ii]; |
---|
142 | if ( a ) delete a; |
---|
143 | } |
---|
144 | partialSumSigma.clear(); |
---|
145 | } |
---|
146 | // fill new data |
---|
147 | if (numOfCouples>0) { |
---|
148 | for (G4int i=0; i<numOfCouples; i++) { |
---|
149 | G4double cute = DBL_MAX; |
---|
150 | |
---|
151 | // protection for usage with extrapolator |
---|
152 | if(i < nc) cute = cuts[i]; |
---|
153 | |
---|
154 | const G4MaterialCutsCouple* couple = |
---|
155 | theCoupleTable->GetMaterialCutsCouple(i); |
---|
156 | const G4Material* material = couple->GetMaterial(); |
---|
157 | G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute); |
---|
158 | partialSumSigma.push_back(dv); |
---|
159 | } |
---|
160 | } |
---|
161 | } |
---|
162 | |
---|
163 | // define pointer to G4ParticleChange |
---|
164 | if(!fParticleChange) fParticleChange = GetParticleChangeForLoss(); |
---|
165 | } |
---|
166 | |
---|
167 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
168 | |
---|
169 | G4double G4MuBremsstrahlungModel::ComputeDEDXPerVolume( |
---|
170 | const G4Material* material, |
---|
171 | const G4ParticleDefinition*, |
---|
172 | G4double kineticEnergy, |
---|
173 | G4double cutEnergy) |
---|
174 | { |
---|
175 | G4double dedx = 0.0; |
---|
176 | if (kineticEnergy <= lowestKinEnergy) return dedx; |
---|
177 | |
---|
178 | G4double tmax = kineticEnergy; |
---|
179 | G4double cut = std::min(cutEnergy,tmax); |
---|
180 | if(cut < minThreshold) cut = minThreshold; |
---|
181 | |
---|
182 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
183 | const G4double* theAtomicNumDensityVector = |
---|
184 | material->GetAtomicNumDensityVector(); |
---|
185 | |
---|
186 | // loop for elements in the material |
---|
187 | for (size_t i=0; i<material->GetNumberOfElements(); i++) { |
---|
188 | |
---|
189 | G4double loss = |
---|
190 | ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut); |
---|
191 | |
---|
192 | dedx += loss*theAtomicNumDensityVector[i]; |
---|
193 | } |
---|
194 | // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl; |
---|
195 | if(dedx < 0.) dedx = 0.; |
---|
196 | return dedx; |
---|
197 | } |
---|
198 | |
---|
199 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
200 | |
---|
201 | G4double G4MuBremsstrahlungModel::ComputMuBremLoss(G4double Z, |
---|
202 | G4double tkin, G4double cut) |
---|
203 | { |
---|
204 | G4double totalEnergy = mass + tkin; |
---|
205 | G4double ak1 = 0.05; |
---|
206 | G4int k2=5; |
---|
207 | G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623}; |
---|
208 | G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566}; |
---|
209 | G4double loss = 0.; |
---|
210 | |
---|
211 | G4double vcut = cut/totalEnergy; |
---|
212 | G4double vmax = tkin/totalEnergy; |
---|
213 | |
---|
214 | G4double aaa = 0.; |
---|
215 | G4double bbb = vcut; |
---|
216 | if(vcut>vmax) bbb=vmax ; |
---|
217 | G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ; |
---|
218 | G4double hhh=(bbb-aaa)/float(kkk) ; |
---|
219 | |
---|
220 | G4double aa = aaa; |
---|
221 | for(G4int l=0; l<kkk; l++) |
---|
222 | { |
---|
223 | for(G4int i=0; i<6; i++) |
---|
224 | { |
---|
225 | G4double ep = (aa + xgi[i]*hhh)*totalEnergy; |
---|
226 | loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep); |
---|
227 | } |
---|
228 | aa += hhh; |
---|
229 | } |
---|
230 | |
---|
231 | loss *=hhh*totalEnergy ; |
---|
232 | |
---|
233 | return loss; |
---|
234 | } |
---|
235 | |
---|
236 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
237 | |
---|
238 | G4double G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection( |
---|
239 | G4double tkin, |
---|
240 | G4double Z, |
---|
241 | G4double cut) |
---|
242 | { |
---|
243 | G4double totalEnergy = tkin + mass; |
---|
244 | G4double ak1 = 2.3; |
---|
245 | G4int k2 = 4; |
---|
246 | G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623}; |
---|
247 | G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566}; |
---|
248 | G4double cross = 0.; |
---|
249 | |
---|
250 | if(cut >= tkin) return cross; |
---|
251 | |
---|
252 | G4double vcut = cut/totalEnergy; |
---|
253 | G4double vmax = tkin/totalEnergy; |
---|
254 | |
---|
255 | G4double aaa = log(vcut); |
---|
256 | G4double bbb = log(vmax); |
---|
257 | G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ; |
---|
258 | G4double hhh = (bbb-aaa)/G4double(kkk); |
---|
259 | |
---|
260 | G4double aa = aaa; |
---|
261 | |
---|
262 | for(G4int l=0; l<kkk; l++) |
---|
263 | { |
---|
264 | for(G4int i=0; i<6; i++) |
---|
265 | { |
---|
266 | G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy; |
---|
267 | cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep); |
---|
268 | } |
---|
269 | aa += hhh; |
---|
270 | } |
---|
271 | |
---|
272 | cross *=hhh; |
---|
273 | |
---|
274 | //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl; |
---|
275 | |
---|
276 | return cross; |
---|
277 | } |
---|
278 | |
---|
279 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
280 | |
---|
281 | G4double G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection( |
---|
282 | G4double tkin, |
---|
283 | G4double Z, |
---|
284 | G4double gammaEnergy) |
---|
285 | // differential cross section |
---|
286 | { |
---|
287 | G4double dxsection = 0.; |
---|
288 | |
---|
289 | if( gammaEnergy > tkin) return dxsection ; |
---|
290 | |
---|
291 | G4double E = tkin + mass ; |
---|
292 | G4double v = gammaEnergy/E ; |
---|
293 | G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ; |
---|
294 | G4double rab0=delta*sqrte ; |
---|
295 | |
---|
296 | G4int iz = G4int(Z); |
---|
297 | if(iz < 1) iz = 1; |
---|
298 | |
---|
299 | G4double z13 = 1.0/nist->GetZ13(iz); |
---|
300 | G4double dn = 1.54*nist->GetA27(iz); |
---|
301 | |
---|
302 | G4double b,b1,dnstar ; |
---|
303 | |
---|
304 | if(1 == iz) |
---|
305 | { |
---|
306 | b = bh; |
---|
307 | b1 = bh1; |
---|
308 | dnstar = dn; |
---|
309 | } |
---|
310 | else |
---|
311 | { |
---|
312 | b = btf; |
---|
313 | b1 = btf1; |
---|
314 | dnstar = dn/std::pow(dn, 1./Z); |
---|
315 | } |
---|
316 | |
---|
317 | // nucleus contribution logarithm |
---|
318 | G4double rab1=b*z13; |
---|
319 | G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))* |
---|
320 | (mass+delta*(dnstar*sqrte-2.))) ; |
---|
321 | if(fn <0.) fn = 0. ; |
---|
322 | // electron contribution logarithm |
---|
323 | G4double epmax1=E/(1.+0.5*mass*rmass/E) ; |
---|
324 | G4double fe=0.; |
---|
325 | if(gammaEnergy<epmax1) |
---|
326 | { |
---|
327 | G4double rab2=b1*z13*z13 ; |
---|
328 | fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))* |
---|
329 | (electron_mass_c2+rab0*rab2))) ; |
---|
330 | if(fe<0.) fe=0. ; |
---|
331 | } |
---|
332 | |
---|
333 | dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy; |
---|
334 | |
---|
335 | return dxsection; |
---|
336 | } |
---|
337 | |
---|
338 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
339 | |
---|
340 | G4double G4MuBremsstrahlungModel::ComputeCrossSectionPerAtom( |
---|
341 | const G4ParticleDefinition*, |
---|
342 | G4double kineticEnergy, |
---|
343 | G4double Z, G4double, |
---|
344 | G4double cutEnergy, |
---|
345 | G4double maxEnergy) |
---|
346 | { |
---|
347 | G4double cross = 0.0; |
---|
348 | if (kineticEnergy <= lowestKinEnergy) return cross; |
---|
349 | G4double tmax = std::min(maxEnergy, kineticEnergy); |
---|
350 | G4double cut = std::min(cutEnergy, kineticEnergy); |
---|
351 | if(cut < minThreshold) cut = minThreshold; |
---|
352 | if (cut >= tmax) return cross; |
---|
353 | |
---|
354 | cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut); |
---|
355 | if(tmax < kineticEnergy) { |
---|
356 | cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax); |
---|
357 | } |
---|
358 | return cross; |
---|
359 | } |
---|
360 | |
---|
361 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
362 | |
---|
363 | G4DataVector* G4MuBremsstrahlungModel::ComputePartialSumSigma( |
---|
364 | const G4Material* material, |
---|
365 | G4double kineticEnergy, |
---|
366 | G4double cut) |
---|
367 | |
---|
368 | // Build the table of cross section per element. |
---|
369 | // The table is built for material |
---|
370 | // This table is used to select randomly an element in the material. |
---|
371 | { |
---|
372 | G4int nElements = material->GetNumberOfElements(); |
---|
373 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
374 | const G4double* theAtomNumDensityVector = |
---|
375 | material->GetAtomicNumDensityVector(); |
---|
376 | |
---|
377 | G4DataVector* dv = new G4DataVector(); |
---|
378 | |
---|
379 | G4double cross = 0.0; |
---|
380 | |
---|
381 | for (G4int i=0; i<nElements; i++ ) { |
---|
382 | cross += theAtomNumDensityVector[i] |
---|
383 | * ComputeMicroscopicCrossSection(kineticEnergy, |
---|
384 | (*theElementVector)[i]->GetZ(), cut); |
---|
385 | dv->push_back(cross); |
---|
386 | } |
---|
387 | return dv; |
---|
388 | } |
---|
389 | |
---|
390 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
391 | |
---|
392 | void G4MuBremsstrahlungModel::SampleSecondaries( |
---|
393 | std::vector<G4DynamicParticle*>* vdp, |
---|
394 | const G4MaterialCutsCouple* couple, |
---|
395 | const G4DynamicParticle* dp, |
---|
396 | G4double minEnergy, |
---|
397 | G4double maxEnergy) |
---|
398 | { |
---|
399 | G4double kineticEnergy = dp->GetKineticEnergy(); |
---|
400 | // check against insufficient energy |
---|
401 | G4double tmax = std::min(kineticEnergy, maxEnergy); |
---|
402 | G4double tmin = std::min(kineticEnergy, minEnergy); |
---|
403 | if(tmin < minThreshold) tmin = minThreshold; |
---|
404 | if(tmin >= tmax) return; |
---|
405 | |
---|
406 | // ===== sampling of energy transfer ====== |
---|
407 | |
---|
408 | G4ParticleMomentum partDirection = dp->GetMomentumDirection(); |
---|
409 | |
---|
410 | // select randomly one element constituing the material |
---|
411 | const G4Element* anElement = SelectRandomAtom(couple); |
---|
412 | G4double Z = anElement->GetZ(); |
---|
413 | |
---|
414 | G4double totalEnergy = kineticEnergy + mass; |
---|
415 | G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass)); |
---|
416 | |
---|
417 | G4double func1 = tmin* |
---|
418 | ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin); |
---|
419 | |
---|
420 | G4double lnepksi, epksi; |
---|
421 | G4double func2; |
---|
422 | |
---|
423 | do { |
---|
424 | lnepksi = log(tmin) + G4UniformRand()*log(kineticEnergy/tmin); |
---|
425 | epksi = exp(lnepksi); |
---|
426 | func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi); |
---|
427 | |
---|
428 | } while(func2 < func1*G4UniformRand()); |
---|
429 | |
---|
430 | G4double gEnergy = epksi; |
---|
431 | |
---|
432 | // ===== sample angle ===== |
---|
433 | |
---|
434 | G4double gam = totalEnergy/mass; |
---|
435 | G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0); |
---|
436 | G4double rmax2= rmax*rmax; |
---|
437 | G4double x = G4UniformRand()*rmax2/(1.0 + rmax2); |
---|
438 | |
---|
439 | G4double theta = sqrt(x/(1.0 - x))/gam; |
---|
440 | G4double sint = sin(theta); |
---|
441 | G4double phi = twopi * G4UniformRand() ; |
---|
442 | G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ; |
---|
443 | |
---|
444 | G4ThreeVector gDirection(dirx, diry, dirz); |
---|
445 | gDirection.rotateUz(partDirection); |
---|
446 | |
---|
447 | partDirection *= totalMomentum; |
---|
448 | partDirection -= gEnergy*gDirection; |
---|
449 | partDirection = partDirection.unit(); |
---|
450 | |
---|
451 | // primary change |
---|
452 | kineticEnergy -= gEnergy; |
---|
453 | fParticleChange->SetProposedKineticEnergy(kineticEnergy); |
---|
454 | fParticleChange->SetProposedMomentumDirection(partDirection); |
---|
455 | |
---|
456 | // save secondary |
---|
457 | G4DynamicParticle* aGamma = |
---|
458 | new G4DynamicParticle(theGamma,gDirection,gEnergy); |
---|
459 | vdp->push_back(aGamma); |
---|
460 | } |
---|
461 | |
---|
462 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
463 | |
---|
464 | const G4Element* G4MuBremsstrahlungModel::SelectRandomAtom( |
---|
465 | const G4MaterialCutsCouple* couple) const |
---|
466 | { |
---|
467 | // select randomly 1 element within the material |
---|
468 | |
---|
469 | const G4Material* material = couple->GetMaterial(); |
---|
470 | G4int nElements = material->GetNumberOfElements(); |
---|
471 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
472 | if(1 == nElements) return (*theElementVector)[0]; |
---|
473 | else if(1 > nElements) return 0; |
---|
474 | |
---|
475 | G4DataVector* dv = partialSumSigma[couple->GetIndex()]; |
---|
476 | G4double rval = G4UniformRand()*((*dv)[nElements-1]); |
---|
477 | for (G4int i=0; i<nElements; i++) { |
---|
478 | if (rval <= (*dv)[i]) return (*theElementVector)[i]; |
---|
479 | } |
---|
480 | return (*theElementVector)[nElements-1]; |
---|
481 | } |
---|
482 | |
---|
483 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|