Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4WentzelVIModel.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 // $Id: G4WentzelVIModel.cc 66596 2012-12-23 14:57:45Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4WentzelVIModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 09.04.2008 from G4MuMscModel
38 //
39 // Modifications:
40 // 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
41 // compute cross sections and sample scattering angle
42 //
43 //
44 // Class Description:
45 //
46 // Implementation of the model of multiple scattering based on
47 // G.Wentzel, Z. Phys. 40 (1927) 590.
48 // H.W.Lewis, Phys Rev 78 (1950) 526.
49 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
50 // L.Urban, CERN-OPEN-2006-077.
51 
52 // -------------------------------------------------------------------
53 //
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
58 #include "G4WentzelVIModel.hh"
59 #include "G4PhysicalConstants.hh"
60 #include "G4SystemOfUnits.hh"
61 #include "Randomize.hh"
63 #include "G4PhysicsTableHelper.hh"
64 #include "G4ElementVector.hh"
65 #include "G4ProductionCutsTable.hh"
66 #include "G4LossTableManager.hh"
67 #include "G4Pow.hh"
68 #include "G4NistManager.hh"
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
72 using namespace std;
73 
75  G4VMscModel(nam),
76  numlimit(0.1),
77  currentCouple(0),
78  cosThetaMin(1.0),
79  inside(false),
80  singleScatteringMode(false)
81 {
82  invsqrt12 = 1./sqrt(12.);
83  tlimitminfix = 1.e-6*mm;
84  lowEnergyLimit = 1.0*eV;
85  particle = 0;
86  nelments = 5;
87  xsecn.resize(nelments);
88  prob.resize(nelments);
89  theManager = G4LossTableManager::Instance();
90  fNistManager = G4NistManager::Instance();
91  fG4pow = G4Pow::GetInstance();
92  wokvi = new G4WentzelOKandVIxSection();
93 
94  preKinEnergy = tPathLength = zPathLength = lambdaeff = currentRange = xtsec = 0;
95  currentMaterialIndex = 0;
96  cosThetaMax = cosTetMaxNuc = 1.0;
97 
98  fParticleChange = 0;
99  currentCuts = 0;
100  currentMaterial = 0;
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
106 {
107  delete wokvi;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111 
113  const G4DataVector& cuts)
114 {
115  // reset parameters
116  SetupParticle(p);
117  currentRange = 0.0;
118 
119  cosThetaMax = cos(PolarAngleLimit());
120  wokvi->Initialise(p, cosThetaMax);
121  /*
122  G4cout << "G4WentzelVIModel: " << particle->GetParticleName()
123  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
124  << G4endl;
125  */
126  currentCuts = &cuts;
127 
128  // set values of some data members
129  fParticleChange = GetParticleChangeForMSC(p);
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
135  const G4ParticleDefinition* p,
136  G4double kinEnergy,
138  G4double cutEnergy, G4double)
139 {
140  G4double cross = 0.0;
141  if(p != particle) { SetupParticle(p); }
142  if(kinEnergy < lowEnergyLimit) { return cross; }
143  if(!CurrentCouple()) {
144  G4Exception("G4WentzelVIModel::ComputeCrossSectionPerAtom", "em0011",
145  FatalException, " G4MaterialCutsCouple is not defined");
146  return 0.0;
147  }
148  DefineMaterial(CurrentCouple());
149  cosTetMaxNuc = wokvi->SetupKinematic(kinEnergy, currentMaterial);
150  if(cosTetMaxNuc < 1.0) {
151  cosTetMaxNuc = wokvi->SetupTarget(G4lrint(Z), cutEnergy);
152  cross = wokvi->ComputeTransportCrossSectionPerAtom(cosTetMaxNuc);
153  /*
154  if(p->GetParticleName() == "e-")
155  G4cout << "G4WentzelVIModel::CS: Z= " << G4int(Z) << " e(MeV)= " << kinEnergy
156  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
157  << " " << particle->GetParticleName() << G4endl;
158  */
159  }
160  return cross;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164 
166 {
167  SetupParticle(track->GetDynamicParticle()->GetDefinition());
168  inside = false;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172 
174  const G4Track& track,
175  G4double& currentMinimalStep)
176 {
177  G4double tlimit = currentMinimalStep;
178  const G4DynamicParticle* dp = track.GetDynamicParticle();
179  G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
180  G4StepStatus stepStatus = sp->GetStepStatus();
181  singleScatteringMode = false;
182  //G4cout << "G4WentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
183  // << stepStatus << G4endl;
184 
185  // initialisation for each step, lambda may be computed from scratch
186  preKinEnergy = dp->GetKineticEnergy();
187  DefineMaterial(track.GetMaterialCutsCouple());
188  lambdaeff = GetTransportMeanFreePath(particle,preKinEnergy);
189  currentRange = GetRange(particle,preKinEnergy,currentCouple);
190  cosTetMaxNuc = wokvi->SetupKinematic(preKinEnergy, currentMaterial);
191 
192  //G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
193  // << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
194 
195  // extra check for abnormal situation
196  // this check needed to run MSC with eIoni and eBrem inactivated
197  if(tlimit > currentRange) { tlimit = currentRange; }
198 
199  // stop here if small range particle
200  if(inside || tlimit < tlimitminfix) {
201  return ConvertTrueToGeom(tlimit, currentMinimalStep);
202  }
203 
204  // pre step
205  G4double presafety = sp->GetSafety();
206  // far from geometry boundary
207  if(currentRange < presafety) {
208  inside = true;
209  return ConvertTrueToGeom(tlimit, currentMinimalStep);
210  }
211 
212  // compute presafety again if presafety <= 0 and no boundary
213  // i.e. when it is needed for optimization purposes
214  if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
215  presafety = ComputeSafety(sp->GetPosition(), tlimit);
216  if(currentRange < presafety) {
217  inside = true;
218  return ConvertTrueToGeom(tlimit, currentMinimalStep);
219  }
220  }
221  /*
222  G4cout << "e(MeV)= " << preKinEnergy/MeV
223  << " " << particle->GetParticleName()
224  << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
225  << " R(mm)= " <<currentRange/mm
226  << " L0(mm^-1)= " << lambdaeff*mm
227  <<G4endl;
228  */
229 
230  // natural limit for high energy
231  G4double rlimit = std::max(facrange*currentRange,
232  0.7*(1.0 - cosTetMaxNuc)*lambdaeff);
233 
234  // low-energy e-
235  if(cosThetaMax > cosTetMaxNuc) {
236  rlimit = std::min(rlimit, facsafety*presafety);
237  }
238 
239  // cut correction
240  G4double rcut = currentCouple->GetProductionCuts()->GetProductionCut(1);
241  //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
242  // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
243  //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
244  if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
245 
246  if(rlimit < tlimit) { tlimit = rlimit; }
247 
248  tlimit = std::max(tlimit, tlimitminfix);
249 
250  // step limit in infinite media
251  tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
252 
253  //compute geomlimit and force few steps within a volume
255  && stepStatus == fGeomBoundary) {
256 
257  G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
258  tlimit = std::min(tlimit, geomlimit/facgeom);
259  }
260 
261  /*
262  G4cout << particle->GetParticleName() << " e= " << preKinEnergy
263  << " L0= " << lambdaeff << " R= " << currentRange
264  << "tlimit= " << tlimit
265  << " currentMinimalStep= " << currentMinimalStep << G4endl;
266  */
267  return ConvertTrueToGeom(tlimit, currentMinimalStep);
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271 
273 {
274  tPathLength = truelength;
275  zPathLength = tPathLength;
276 
277  if(lambdaeff > 0.0 && lambdaeff != DBL_MAX) {
278  G4double tau = tPathLength/lambdaeff;
279  //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
280  // << " Leff= " << lambdaeff << " tau= " << tau << G4endl;
281  // small step
282  if(tau < numlimit) {
283  zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
284 
285  // medium step
286  } else {
287  G4double e1 = 0.0;
288  if(currentRange > tPathLength) {
289  e1 = GetEnergy(particle,currentRange-tPathLength,currentCouple);
290  }
291  e1 = 0.5*(e1 + preKinEnergy);
292  cosTetMaxNuc = wokvi->SetupKinematic(e1, currentMaterial);
293  lambdaeff = GetTransportMeanFreePath(particle,e1);
294  zPathLength = lambdaeff*(1.0 - exp(-tPathLength/lambdaeff));
295  }
296  } else { lambdaeff = DBL_MAX; }
297  //G4cout<<"Comp.geom: zLength= "<<zPathLength<<" tLength= "<<tPathLength<<G4endl;
298  return zPathLength;
299 }
300 
301 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
302 
304 {
305  // initialisation of single scattering x-section
306  xtsec = 0.0;
307  cosThetaMin = cosTetMaxNuc;
308  /*
309  G4cout << "ComputeTrueStepLength: Step= " << geomStepLength
310  << " Lambda= " << lambdaeff
311  << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
312  */
313  // pathalogical case
314  if(lambdaeff == DBL_MAX) {
315  singleScatteringMode = true;
316  zPathLength = geomStepLength;
317  tPathLength = geomStepLength;
318  cosThetaMin = 1.0;
319 
320  // normal case
321  } else {
322 
323  // small step use only single scattering
324  const G4double singleScatLimit = 1.0e-7;
325  if(geomStepLength < lambdaeff*singleScatLimit*(1.0 - cosTetMaxNuc)) {
326  singleScatteringMode = true;
327  cosThetaMin = 1.0;
328  zPathLength = geomStepLength;
329  tPathLength = geomStepLength;
330 
331  // step defined by transportation
332  } else if(geomStepLength != zPathLength) {
333 
334  // step defined by transportation
335  zPathLength = geomStepLength;
336  G4double tau = geomStepLength/lambdaeff;
337  tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
338 
339  // energy correction for a big step
340  if(tau > numlimit) {
341  G4double e1 = 0.0;
342  if(currentRange > tPathLength) {
343  e1 = GetEnergy(particle,currentRange-tPathLength,currentCouple);
344  }
345  e1 = 0.5*(e1 + preKinEnergy);
346  cosTetMaxNuc = wokvi->SetupKinematic(e1, currentMaterial);
347  lambdaeff = GetTransportMeanFreePath(particle,e1);
348  tau = zPathLength/lambdaeff;
349 
350  if(tau < 0.999999) { tPathLength = -lambdaeff*log(1.0 - tau); }
351  else { tPathLength = currentRange; }
352  }
353  }
354  }
355 
356  // check of step length
357  // define threshold angle between single and multiple scattering
358  if(!singleScatteringMode) { cosThetaMin = 1.0 - 1.5*tPathLength/lambdaeff; }
359 
360  // recompute transport cross section - do not change energy
361  // anymore - cannot be applied for big steps
362  if(cosThetaMin > cosTetMaxNuc) {
363 
364  // new computation
365  G4double cross = ComputeXSectionPerVolume();
366  //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec << G4endl;
367  if(cross <= 0.0) {
368  singleScatteringMode = true;
369  tPathLength = zPathLength;
370  lambdaeff = DBL_MAX;
371  cosThetaMin = 1.0;
372  } else if(xtsec > 0.0) {
373 
374  lambdaeff = 1./cross;
375  G4double tau = zPathLength*cross;
376  if(tau < numlimit) { tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0); }
377  else if(tau < 0.999999) { tPathLength = -lambdaeff*log(1.0 - tau); }
378  else { tPathLength = currentRange; }
379 
380  if(tPathLength > currentRange) { tPathLength = currentRange; }
381  }
382  }
383 
384  /*
385  G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
386  <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
387  G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
388  << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
389  << " e(MeV)= " << preKinEnergy/MeV << " "
390  << singleScatteringMode << G4endl;
391  */
392  return tPathLength;
393 }
394 
395 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
396 
399  G4double safety)
400 {
401  fDisplacement.set(0.0,0.0,0.0);
402  //G4cout << "!##! G4WentzelVIModel::SampleScattering for "
403  // << particle->GetParticleName() << G4endl;
404 
405  // ignore scattering for zero step length and energy below the limit
406  if(preKinEnergy < lowEnergyLimit || tPathLength <= 0.0)
407  { return fDisplacement; }
408 
409  G4double invlambda = 0.0;
410  if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
411 
412  // use average kinetic energy over the step
413  G4double cut = (*currentCuts)[currentMaterialIndex];
414  /*
415  G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
416  << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
417  << " xmsc= " << tPathLength*invlambda
418  << " safety= " << safety << G4endl;
419  */
420 
421  // step limit due msc
422  G4double x0 = tPathLength;
423  // const G4double thinlimit = 0.5;
424  const G4double thinlimit = 0.1;
425  G4int nMscSteps = 1;
426  // large scattering angle case - two step approach
427  if(!singleScatteringMode && tPathLength*invlambda > thinlimit
428  && safety > tlimitminfix) {
429  x0 *= 0.5;
430  nMscSteps = 2;
431  }
432 
433  // step limit due to single scattering
434  G4double x1 = 2*tPathLength;
435  if(0.0 < xtsec) { x1 = -log(G4UniformRand())/xtsec; }
436 
437  // no scattering case
438  if(singleScatteringMode && x1 > tPathLength)
439  { return fDisplacement; }
440 
441  const G4ElementVector* theElementVector =
442  currentMaterial->GetElementVector();
443  G4int nelm = currentMaterial->GetNumberOfElements();
444 
445  // geometry
446  G4double sint, cost, phi;
447  G4ThreeVector temp(0.0,0.0,1.0);
448 
449  // current position and direction relative to the end point
450  // because of magnetic field geometry is computed relatively to the
451  // end point of the step
452  G4ThreeVector dir(0.0,0.0,1.0);
453  fDisplacement.set(0.0,0.0,-zPathLength);
454 
455  G4double mscfac = zPathLength/tPathLength;
456 
457  // start a loop
458  G4double x2 = x0;
459  G4double step, z;
460  G4bool singleScat;
461  /*
462  G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
463  << " 1-cost1= " << 1 - cosThetaMin << " " << singleScatteringMode
464  << " xtsec= " << xtsec << G4endl;
465  */
466  do {
467 
468  // single scattering case
469  if(singleScatteringMode || x1 <= x2) {
470  step = x1;
471  singleScat = true;
472  } else {
473  step = x2;
474  singleScat = false;
475  }
476 
477  // new position
478  fDisplacement += step*mscfac*dir;
479 
480  if(singleScat) {
481 
482  // select element
483  G4int i = 0;
484  if(nelm > 1) {
485  G4double qsec = G4UniformRand()*xtsec;
486  for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
487  }
488  G4double cosTetM =
489  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
490  //G4cout << "!!! " << cosThetaMin << " " << cosTetM << " " << prob[i] << G4endl;
491  temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
492 
493  // direction is changed
494  temp.rotateUz(dir);
495  dir = temp;
496 
497  // new proposed step length
498  x2 -= step;
499  x1 = -log(G4UniformRand())/xtsec;
500  if(singleScatteringMode && x1 > x2) { break; }
501 
502  // multiple scattering
503  } else {
504  --nMscSteps;
505  x1 -= step;
506  x2 = x0;
507 
508  // for multiple scattering x0 should be used as a step size
509  G4double z0 = x0*invlambda;
510 
511  // sample z in interval 0 - 1
512  do {
513  G4double zzz = 0.0;
514  if(z0 > 0.1) { zzz = exp(-1.0/z0); }
515  z = -z0*log(1.0 - (1.0 - zzz)*G4UniformRand());
516  } while(z > 1.0);
517  cost = 1.0 - 2.0*z/*factCM*/;
518  if(cost > 1.0) { cost = 1.0; }
519  else if(cost < -1.0) { cost =-1.0; }
520  sint = sqrt((1.0 - cost)*(1.0 + cost));
521  phi = twopi*G4UniformRand();
522  G4double vx1 = sint*cos(phi);
523  G4double vy1 = sint*sin(phi);
524 
525  // change direction
526  temp.set(vx1,vy1,cost);
527  temp.rotateUz(dir);
528  dir = temp;
529 
530  // lateral displacement
531  if (latDisplasment && safety > tlimitminfix) {
532  G4double rms = invsqrt12*sqrt(2*z0);
533  G4double r = x0*mscfac;
534  G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(0.0,1.0));
535  G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(0.0,1.0));
536  G4double dz;
537  G4double d = r*r - dx*dx - dy*dy;
538  if(d >= 0.0) { dz = sqrt(d) - r; }
539  else { dx = dy = dz = 0.0; }
540 
541  // change position
542  temp.set(dx,dy,dz);
543  temp.rotateUz(dir);
544  fDisplacement += temp;
545  }
546  }
547  } while (0 < nMscSteps);
548 
549  dir.rotateUz(oldDirection);
550 
551  //G4cout << "G4WentzelVIModel sampling of scattering is done" << G4endl;
552  // end of sampling -------------------------------
553 
554  fParticleChange->ProposeMomentumDirection(dir);
555 
556  // lateral displacement
557  fDisplacement.rotateUz(oldDirection);
558 
559  /*
560  G4cout << " r(mm)= " << fDisplacement.mag()
561  << " safety= " << safety
562  << " trueStep(mm)= " << tPathLength
563  << " geomStep(mm)= " << zPathLength
564  << " x= " << fDisplacement.x()
565  << " y= " << fDisplacement.y()
566  << " z= " << fDisplacement.z()
567  << G4endl;
568  */
569 
570  //G4cout<< "G4WentzelVIModel::SampleScattering end NewDir= " << dir<< G4endl;
571  return fDisplacement;
572 }
573 
574 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
575 
576 G4double G4WentzelVIModel::ComputeXSectionPerVolume()
577 {
578  // prepare recomputation of x-sections
579  const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
580  const G4double* theAtomNumDensityVector =
581  currentMaterial->GetVecNbOfAtomsPerVolume();
582  G4int nelm = currentMaterial->GetNumberOfElements();
583  if(nelm > nelments) {
584  nelments = nelm;
585  xsecn.resize(nelm);
586  prob.resize(nelm);
587  }
588  G4double cut = (*currentCuts)[currentMaterialIndex];
589  // cosTetMaxNuc = wokvi->GetCosThetaNuc();
590 
591  // check consistency
592  xtsec = 0.0;
593  if(cosTetMaxNuc > cosThetaMin) { return 0.0; }
594 
595  // loop over elements
596  G4double xs = 0.0;
597  for (G4int i=0; i<nelm; ++i) {
598  G4double costm =
599  wokvi->SetupTarget(G4lrint((*theElementVector)[i]->GetZ()), cut);
600  G4double density = theAtomNumDensityVector[i];
601 
602  G4double esec = 0.0;
603  if(costm < cosThetaMin) {
604 
605  // recompute the transport x-section
606  if(1.0 > cosThetaMin) {
607  xs += density*wokvi->ComputeTransportCrossSectionPerAtom(cosThetaMin);
608  }
609  // recompute the total x-section
610  G4double nucsec = wokvi->ComputeNuclearCrossSection(cosThetaMin, costm);
611  esec = wokvi->ComputeElectronCrossSection(cosThetaMin, costm);
612  nucsec += esec;
613  if(nucsec > 0.0) { esec /= nucsec; }
614  xtsec += nucsec*density;
615  }
616  xsecn[i] = xtsec;
617  prob[i] = esec;
618  //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
619  // << " 1-cosThetaMin= " << 1-cosThetaMin
620  // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
621  }
622 
623  //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
624  // << " txsec(1/mm)= " << xtsec <<G4endl;
625  return xs;
626 }
627 
628 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......