Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4WentzelVIModel
33//
34// Author: V.Ivanchenko
35//
36// Creation date: 09.04.2008 from G4MuMscModel
37//
38// Modifications:
39// 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
40// compute cross sections and sample scattering angle
41//
42//
43// Class Description:
44//
45// Implementation of the model of multiple scattering based on
46// G.Wentzel, Z. Phys. 40 (1927) 590.
47// H.W.Lewis, Phys Rev 78 (1950) 526.
48// J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
49// L.Urban, CERN-OPEN-2006-077.
50
51// -------------------------------------------------------------------
52//
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
57#include "G4WentzelVIModel.hh"
59#include "G4SystemOfUnits.hh"
60#include "Randomize.hh"
63#include "G4ElementVector.hh"
65#include "G4EmParameters.hh"
66#include "G4Log.hh"
67#include "G4Exp.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71using namespace std;
72
73const G4double invsqrt12 = 1./std::sqrt(12.);
74const G4double numlimit = 0.1;
76
78 : G4VMscModel(nam),
79 singleScatteringMode(false),
80 isCombined(comb),
81 useSecondMoment(false)
82{
83 tlimitminfix = 1.e-6*CLHEP::mm;
84 lowEnergyLimit = 1.0*CLHEP::eV;
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
92{
93 delete wokvi;
94 if(IsMaster()) {
95 delete fSecondMoments;
96 fSecondMoments = nullptr;
97 }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
103 const G4DataVector& cuts)
104{
105 // reset parameters
106 SetupParticle(p);
108 currentRange = 0.0;
109
110 if(isCombined) {
112 if(tet <= 0.0) { cosThetaMax = 1.0; }
113 else if(tet < CLHEP::pi) { cosThetaMax = cos(tet); }
114 }
115 //G4cout << "G4WentzelVIModel::Initialise " << p->GetParticleName()
116 // << " " << this << " " << wokvi << G4endl;
117
119 /*
120 G4cout << "G4WentzelVIModel: " << particle->GetParticleName()
121 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMax
122 << " SingScatFactor= " << ssFactor
123 << G4endl;
124 */
125 currentCuts = &cuts;
126
127 // set values of some data members
129
130 // Access to materials
131 const G4ProductionCutsTable* theCoupleTable =
133 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
134 nelments = 0;
135 for(G4int i=0; i<numOfCouples; ++i) {
136 G4int nelm = (G4int)theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetNumberOfElements();
137 nelments = std::max(nelments, nelm);
138 }
139 xsecn.resize(nelments);
140 prob.resize(nelments);
141
142 // build second moment table only if transport table is build
144 if(useSecondMoment && IsMaster() && nullptr != table) {
145
146 //G4cout << "### G4WentzelVIModel::Initialise: build 2nd moment table "
147 // << table << G4endl;
150
151 G4bool splineFlag = true;
152 G4PhysicsVector* aVector = nullptr;
153 G4PhysicsVector* bVector = nullptr;
156 if(emin < emax) {
158 *G4lrint(std::log10(emax/emin));
159 if(n < 3) { n = 3; }
160
161 for(G4int i=0; i<numOfCouples; ++i) {
162
163 //G4cout<< "i= " << i << " Flag= " << fSecondMoments->GetFlag(i)
164 // << G4endl;
165 if(fSecondMoments->GetFlag(i)) {
166 DefineMaterial(theCoupleTable->GetMaterialCutsCouple(i));
167
168 delete (*fSecondMoments)[i];
169 if(nullptr == aVector) {
170 aVector = new G4PhysicsLogVector(emin, emax, n, splineFlag);
171 bVector = aVector;
172 } else {
173 bVector = new G4PhysicsVector(*aVector);
174 }
175 for(std::size_t j=0; j<n; ++j) {
176 G4double e = bVector->Energy(j);
177 bVector->PutValue(j, ComputeSecondMoment(p, e)*e*e);
178 }
179 if(splineFlag) { bVector->FillSecondDerivatives(); }
180 (*fSecondMoments)[i] = bVector;
181 }
182 }
183 }
184 //G4cout << *fSecondMoments << G4endl;
185 }
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189
191 G4VEmModel* masterModel)
192{
193 fSecondMoments = static_cast<G4WentzelVIModel*>(masterModel)
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
212 const G4ParticleDefinition* p,
213 G4double kinEnergy,
215 G4double cutEnergy, G4double)
216{
217 G4double cross = 0.0;
218 SetupParticle(p);
219 if(kinEnergy < lowEnergyLimit) { return cross; }
220 if(nullptr == CurrentCouple()) {
221 G4Exception("G4WentzelVIModel::ComputeCrossSectionPerAtom", "em0011",
222 FatalException, " G4MaterialCutsCouple is not defined");
223 return 0.0;
224 }
227 if(cosTetMaxNuc < 1.0) {
228 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
229 G4double cost = wokvi->SetupTarget(G4lrint(Z), cut);
231 /*
232 if(p->GetParticleName() == "e-")
233 G4cout << "G4WentzelVIModel::CS: Z= " << G4int(Z) << " e(MeV)= "<<kinEnergy
234 << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
235 << " " << particle->GetParticleName() << G4endl;
236 */
237 }
238 return cross;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242
244{
245 /*
246 G4cout << "G4WentzelVIModel::StartTracking " << track << " " << this << " "
247 << track->GetParticleDefinition()->GetParticleName()
248 << " workvi: " << wokvi << G4endl;
249 */
251}
252
253//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
254
256 const G4Track& track,
257 G4double& currentMinimalStep)
258{
259 G4double tlimit = currentMinimalStep;
260 const G4DynamicParticle* dp = track.GetDynamicParticle();
261 const G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
262 G4StepStatus stepStatus = sp->GetStepStatus();
263 singleScatteringMode = false;
264
265 //G4cout << "G4WentzelVIModel::ComputeTruePathLengthLimit stepStatus= "
266 // << stepStatus << " " << track.GetDefinition()->GetParticleName()
267 // << G4endl;
268
269 // initialisation for each step, lambda may be computed from scratch
273 const G4double logPreKinEnergy = dp->GetLogKineticEnergy();
277
278 //G4cout << "lambdaeff= " << lambdaeff << " Range= " << currentRange
279 // << " tlimit= " << tlimit << " 1-cost= " << 1 - cosTetMaxNuc << G4endl;
280
281 // extra check for abnormal situation
282 // this check needed to run MSC with eIoni and eBrem inactivated
283 if(tlimit > currentRange) { tlimit = currentRange; }
284
285 // stop here if small range particle
286 if(tlimit < tlimitminfix) {
287 return ConvertTrueToGeom(tlimit, currentMinimalStep);
288 }
289
290 // pre step
291 G4double presafety = sp->GetSafety();
292 // far from geometry boundary
293 if(currentRange < presafety) {
294 return ConvertTrueToGeom(tlimit, currentMinimalStep);
295 }
296
297 // compute presafety again if presafety <= 0 and no boundary
298 // i.e. when it is needed for optimization purposes
299 if(stepStatus != fGeomBoundary && presafety < tlimitminfix) {
300 presafety = ComputeSafety(sp->GetPosition(), tlimit);
301 if(currentRange < presafety) {
302 return ConvertTrueToGeom(tlimit, currentMinimalStep);
303 }
304 }
305 /*
306 G4cout << "e(MeV)= " << preKinEnergy/MeV
307 << " " << particle->GetParticleName()
308 << " CurLimit(mm)= " << tlimit/mm <<" safety(mm)= " << presafety/mm
309 << " R(mm)= " <<currentRange/mm
310 << " L0(mm^-1)= " << lambdaeff*mm
311 << G4endl;
312 */
313 // natural limit for high energy
314 G4double rlimit = std::max(facrange*currentRange,
316
317 // low-energy e-
319 rlimit = std::min(rlimit, facsafety*presafety);
320 }
321
322 // cut correction
324 //G4cout << "rcut= " << rcut << " rlimit= " << rlimit << " presafety= "
325 // << presafety << " 1-cosThetaMax= " <<1-cosThetaMax
326 //<< " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc << G4endl;
327 if(rcut > rlimit) { rlimit = std::min(rlimit, rcut*sqrt(rlimit/rcut)); }
328
329 tlimit = std::min(tlimit, rlimit);
330 tlimit = std::max(tlimit, tlimitminfix);
331
332 // step limit in infinite media
333 tlimit = std::min(tlimit, 50*currentMaterial->GetRadlen()/facgeom);
334
335 //compute geomlimit and force few steps within a volume
337 && stepStatus == fGeomBoundary) {
338
339 G4double geomlimit = ComputeGeomLimit(track, presafety, currentRange);
340 tlimit = std::min(tlimit, geomlimit/facgeom);
341 }
342 /*
343 G4cout << particle->GetParticleName() << " e= " << preKinEnergy
344 << " L0= " << lambdaeff << " R= " << currentRange
345 << " tlimit= " << tlimit
346 << " currentMinimalStep= " << currentMinimalStep << G4endl;
347 */
348 return ConvertTrueToGeom(tlimit, currentMinimalStep);
349}
350
351//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
352
354{
355 zPathLength = tPathLength = truelength;
356
357 // small step use only single scattering
358 cosThetaMin = 1.0;
360 //G4cout << "xtsec= " << xtsec << " Nav= "
361 // << zPathLength*xtsec << G4endl;
365
366 } else {
367 //G4cout << "ComputeGeomPathLength: tLength= " << tPathLength
368 // << " Leff= " << lambdaeff << G4endl;
369 // small step
372 zPathLength *= (1.0 - 0.5*tau + tau*tau/6.0);
373
374 // medium step
375 } else {
376 G4double e1 = 0.0;
379 }
380 effKinEnergy = 0.5*(e1 + preKinEnergy);
383 //G4cout << " tLength= "<< tPathLength<< " Leff= " << lambdaeff << G4endl;
387 }
388 }
389 }
390 //G4cout << "Comp.geom: zLength= "<<zPathLength<<" tLength= "
391 // << tPathLength<< " Leff= " << lambdaeff << G4endl;
392 return zPathLength;
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
396
398{
399 // initialisation of single scattering x-section
400 /*
401 G4cout << "ComputeTrueStepLength: Step= " << geomStepLength
402 << " geomL= " << zPathLength
403 << " Lambda= " << lambdaeff
404 << " 1-cosThetaMaxNuc= " << 1 - cosTetMaxNuc << G4endl;
405 */
407 zPathLength = tPathLength = geomStepLength;
408
409 } else {
410
411 // step defined by transportation
412 // change both geom and true step lengths
413 if(geomStepLength < zPathLength) {
414
415 // single scattering
416 if(G4int(geomStepLength*xtsec) < minNCollisions) {
417 zPathLength = tPathLength = geomStepLength;
420
421 // multiple scattering
422 } else {
423 // small step
424 if(geomStepLength < numlimit*lambdaeff) {
425 G4double tau = geomStepLength/lambdaeff;
426 tPathLength = geomStepLength*(1.0 + 0.5*tau + tau*tau/3.0);
427
428 // energy correction for a big step
429 } else {
430 tPathLength *= geomStepLength/zPathLength;
431 G4double e1 = 0.0;
434 }
435 effKinEnergy = 0.5*(e1 + preKinEnergy);
438 G4double tau = geomStepLength/lambdaeff;
439
440 if(tau < 0.999999) { tPathLength = -lambdaeff*G4Log(1.0 - tau); }
441 else { tPathLength = currentRange; }
442 }
443 zPathLength = geomStepLength;
444 }
445 }
446 }
447 // check of step length
448 // define threshold angle between single and multiple scattering
451 xtsec = 0.0;
452
453 // recompute transport cross section - do not change energy
454 // anymore - cannot be applied for big steps
456 // new computation
458 //G4cout << "%%%% cross= " << cross << " xtsec= " << xtsec
459 // << " 1-cosTMin= " << 1.0 - cosThetaMin << G4endl;
460 if(cross <= 0.0) {
464 cosThetaMin = 1.0;
465 } else if(xtsec > 0.0) {
466
467 lambdaeff = 1./cross;
468 G4double tau = zPathLength*cross;
469 if(tau < numlimit) {
470 tPathLength = zPathLength*(1.0 + 0.5*tau + tau*tau/3.0);
471 } else if(tau < 0.999999) {
472 tPathLength = -lambdaeff*G4Log(1.0 - tau);
473 } else {
475 }
476 }
477 }
478 }
480 /*
481 G4cout <<"Comp.true: zLength= "<<zPathLength<<" tLength= "<<tPathLength
482 <<" Leff(mm)= "<<lambdaeff/mm<<" sig0(1/mm)= " << xtsec <<G4endl;
483 G4cout << particle->GetParticleName() << " 1-cosThetaMin= " << 1-cosThetaMin
484 << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
485 << " e(MeV)= " << preKinEnergy/MeV << " "
486 << " SSmode= " << singleScatteringMode << G4endl;
487 */
488 return tPathLength;
489}
490
491//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
492
495 G4double /*safety*/)
496{
497 fDisplacement.set(0.0,0.0,0.0);
498 //G4cout << "!##! G4WentzelVIModel::SampleScattering for "
499 // << particle->GetParticleName() << G4endl;
500
501 // ignore scattering for zero step length and energy below the limit
503 { return fDisplacement; }
504
505 G4double invlambda = 0.0;
506 if(lambdaeff < DBL_MAX) { invlambda = 0.5/lambdaeff; }
507
508 // use average kinetic energy over the step
509 G4double cut = (*currentCuts)[currentMaterialIndex];
510 if(fixedCut > 0.0) { cut = fixedCut; }
511 /*
512 G4cout <<"SampleScat: E0(MeV)= "<< preKinEnergy/MeV
513 << " Leff= " << lambdaeff <<" sig0(1/mm)= " << xtsec
514 << " xmsc= " << tPathLength*invlambda
515 << " safety= " << safety << G4endl;
516 */
517 // step limit due msc
518 G4int nMscSteps = 1;
520 G4double z0 = x0*invlambda;
521 //G4double zzz = 0.0;
522 G4double prob2 = 0.0;
523
524 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
525
526 // large scattering angle case - two step approach
528 static const G4double zzmin = 0.05;
529 if(useSecondMoment) {
530 G4double z1 = invlambda*invlambda;
532 prob2 = (z2 - z1)/(1.5*z1 - z2);
533 }
534 // if(z0 > zzmin && safety > tlimitminfix) {
535 if(z0 > zzmin) {
536 x0 *= 0.5;
537 z0 *= 0.5;
538 nMscSteps = 2;
539 }
540 //if(z0 > zzmin) { zzz = G4Exp(-1.0/z0); }
541 G4double zzz = 0.0;
542 if(z0 > zzmin) {
543 zzz = G4Exp(-1.0/z0);
544 z0 += zzz;
545 prob2 *= (1 + zzz);
546 }
547 prob2 /= (1 + prob2);
548 }
549
550 // step limit due to single scattering
551 G4double x1 = 2*tPathLength;
552 if(0.0 < xtsec) { x1 = -G4Log(rndmEngine->flat())/xtsec; }
553
554 // no scattering case
556 { return fDisplacement; }
557
558 const G4ElementVector* theElementVector =
560 std::size_t nelm = currentMaterial->GetNumberOfElements();
561
562 // geometry
563 G4double sint, cost, phi;
564 G4ThreeVector temp(0.0,0.0,1.0);
565
566 // current position and direction relative to the end point
567 // because of magnetic field geometry is computed relatively to the
568 // end point of the step
569 G4ThreeVector dir(0.0,0.0,1.0);
570 fDisplacement.set(0.0,0.0,-zPathLength);
571
573
574 // start a loop
575 G4double x2 = x0;
576 G4double step, z;
577 G4bool singleScat;
578 /*
579 G4cout << "Start of the loop x1(mm)= " << x1 << " x2(mm)= " << x2
580 << " 1-cost1= " << 1 - cosThetaMin << " SSmode= " << singleScatteringMode
581 << " xtsec= " << xtsec << " Nst= " << nMscSteps << G4endl;
582 */
583 do {
584
585 //G4cout << "# x1(mm)= "<< x1<< " x2(mm)= "<< x2 << G4endl;
586 // single scattering case
587 if(singleScatteringMode && x1 > x2) {
588 fDisplacement += x2*mscfac*dir;
589 break;
590 }
591
592 // what is next single of multiple?
593 if(x1 <= x2) {
594 step = x1;
595 singleScat = true;
596 } else {
597 step = x2;
598 singleScat = false;
599 }
600
601 //G4cout << "# step(mm)= "<< step<< " singlScat= "<< singleScat << G4endl;
602
603 // new position
604 fDisplacement += step*mscfac*dir;
605
606 if(singleScat) {
607
608 // select element
609 std::size_t i = 0;
610 if(nelm > 1) {
611 G4double qsec = rndmEngine->flat()*xtsec;
612 for (; i<nelm; ++i) { if(xsecn[i] >= qsec) { break; } }
613 }
614 G4double cosTetM =
615 wokvi->SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
616 //G4cout << "!!! " << cosThetaMin << " " << cosTetM << " "
617 // << prob[i] << G4endl;
618 temp = wokvi->SampleSingleScattering(cosThetaMin, cosTetM, prob[i]);
619
620 // direction is changed
621 temp.rotateUz(dir);
622 dir = temp;
623 //G4cout << dir << G4endl;
624
625 // new proposed step length
626 x2 -= step;
627 x1 = -G4Log(rndmEngine->flat())/xtsec;
628
629 // multiple scattering
630 } else {
631 --nMscSteps;
632 x1 -= step;
633 x2 = x0;
634
635 // sample z in interval 0 - 1
636 G4bool isFirst = true;
637 if(prob2 > 0.0 && rndmEngine->flat() < prob2) { isFirst = false; }
638 do {
639 //z = -z0*G4Log(1.0 - (1.0 - zzz)*rndmEngine->flat());
640 if(isFirst) { z = -G4Log(rndmEngine->flat()); }
641 else { z = G4RandGamma::shoot(rndmEngine, 2.0, 2.0); }
642 z *= z0;
643 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
644 } while(z > 1.0);
645
646 cost = 1.0 - 2.0*z/*factCM*/;
647 if(cost > 1.0) { cost = 1.0; }
648 else if(cost < -1.0) { cost =-1.0; }
649 sint = sqrt((1.0 - cost)*(1.0 + cost));
650 phi = twopi*rndmEngine->flat();
651 G4double vx1 = sint*cos(phi);
652 G4double vy1 = sint*sin(phi);
653
654 // lateral displacement
655 if (latDisplasment) {
656 G4double rms = invsqrt12*sqrt(2*z0);
657 G4double r = x0*mscfac;
658 G4double dx = r*(0.5*vx1 + rms*G4RandGauss::shoot(rndmEngine,0.0,1.0));
659 G4double dy = r*(0.5*vy1 + rms*G4RandGauss::shoot(rndmEngine,0.0,1.0));
660 G4double d = r*r - dx*dx - dy*dy;
661
662 // change position
663 if(d >= 0.0) {
664 temp.set(dx,dy,sqrt(d) - r);
665 temp.rotateUz(dir);
666 fDisplacement += temp;
667 }
668 }
669 // change direction
670 temp.set(vx1,vy1,cost);
671 temp.rotateUz(dir);
672 dir = temp;
673 }
674 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
675 } while (0 < nMscSteps);
676
677 dir.rotateUz(oldDirection);
678
679 //G4cout<<"G4WentzelVIModel sampling is done 1-cost= "<< 1.-dir.z()<<G4endl;
680 // end of sampling -------------------------------
681
683
684 // lateral displacement
685 fDisplacement.rotateUz(oldDirection);
686
687 /*
688 G4cout << " r(mm)= " << fDisplacement.mag()
689 << " safety= " << safety
690 << " trueStep(mm)= " << tPathLength
691 << " geomStep(mm)= " << zPathLength
692 << " x= " << fDisplacement.x()
693 << " y= " << fDisplacement.y()
694 << " z= " << fDisplacement.z()
695 << G4endl;
696 */
697
698 //G4cout<< "G4WentzelVIModel::SampleScattering end NewDir= " << dir<< G4endl;
699 return fDisplacement;
700}
701
702//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
703
705{
706 // prepare recomputation of x-sections
707 const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
708 const G4double* theAtomNumDensityVector =
711 if(nelm > nelments) {
712 nelments = nelm;
713 xsecn.resize(nelm);
714 prob.resize(nelm);
715 }
716
717 // check consistency
718 xtsec = 0.0;
719 if(cosTetMaxNuc >= cosTheta) { return 0.0; }
720
721 G4double cut = (*currentCuts)[currentMaterialIndex];
722 if(fixedCut > 0.0) { cut = fixedCut; }
723
724 // loop over elements
725 G4double xs = 0.0;
726 for (G4int i=0; i<nelm; ++i) {
727 G4double costm =
728 wokvi->SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
729 G4double density = theAtomNumDensityVector[i];
730
731 G4double esec = 0.0;
732 if(costm < cosTheta) {
733
734 // recompute the transport x-section
735 if(1.0 > cosTheta) {
736 xs += density*wokvi->ComputeTransportCrossSectionPerAtom(cosTheta);
737 }
738 // recompute the total x-section
739 G4double nucsec = wokvi->ComputeNuclearCrossSection(cosTheta, costm);
740 esec = wokvi->ComputeElectronCrossSection(cosTheta, costm);
741 nucsec += esec;
742 if(nucsec > 0.0) { esec /= nucsec; }
743 xtsec += nucsec*density;
744 }
745 xsecn[i] = xtsec;
746 prob[i] = esec;
747 //G4cout << i << " xs= " << xs << " xtsec= " << xtsec
748 // << " 1-cosTheta= " << 1-cosTheta
749 // << " 1-cosTetMaxNuc2= " <<1-cosTetMaxNuc2<< G4endl;
750 }
751
752 //G4cout << "ComputeXS result: xsec(1/mm)= " << xs
753 // << " txsec(1/mm)= " << xtsec <<G4endl;
754 return xs;
755}
756
757//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
758
759G4double G4WentzelVIModel::ComputeSecondMoment(const G4ParticleDefinition* p,
760 G4double kinEnergy)
761{
762 G4double xs = 0.0;
763
764 SetupParticle(p);
766
767 if(cosTetMaxNuc >= 1.0) { return xs; }
768
769 const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
770 const G4double* theAtomNumDensityVector =
772 std::size_t nelm = currentMaterial->GetNumberOfElements();
773
774 G4double cut = (*currentCuts)[currentMaterialIndex];
775 if(fixedCut > 0.0) { cut = fixedCut; }
776
777 // loop over elements
778 for (std::size_t i=0; i<nelm; ++i) {
779 G4double costm =
780 wokvi->SetupTarget((*theElementVector)[i]->GetZasInt(), cut);
781 xs += theAtomNumDensityVector[i]
783 }
784 return xs;
785}
786
787//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
788
790{
791 if(val > 0.05) {
792 ssFactor = val;
793 invssFactor = 1.0/(val - 0.05);
794 }
795}
796
797//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
@ fUseDistanceToBoundary
G4StepStatus
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double numlimit
const G4double invsqrt12
const G4double numlimit
const G4int minNCollisions
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
virtual double flat()=0
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4ElementVector * GetElementVector() const
G4double GetRadlen() const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
G4bool GetFlag(std::size_t i) const
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetProductionCut(G4int index) const
G4StepPoint * GetPreStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
G4double PolarAngleLimit() const
G4double LowEnergyLimit() const
G4bool IsMaster() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4double HighEnergyLimit() const
G4double HighEnergyActivationLimit() const
const G4MaterialCutsCouple * CurrentCouple() const
G4PhysicsTable * GetCrossSectionTable()
G4double LowEnergyActivationLimit() const
G4double facrange
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4MscStepLimitType steppingAlgorithm
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
G4bool latDisplasment
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
G4double facsafety
G4ThreeVector fDisplacement
void InitialiseParameters(const G4ParticleDefinition *)
G4double facgeom
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeSecondTransportMoment(G4double CosThetaMax)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
G4double ComputeTransportXSectionPerVolume(G4double cosTheta)
const G4MaterialCutsCouple * currentCouple
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep) override
G4PhysicsTable * GetSecondMomentTable()
std::vector< G4double > prob
G4double SecondMoment(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
~G4WentzelVIModel() override
G4WentzelVIModel(G4bool comb=true, const G4String &nam="WentzelVIUni")
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ComputeTrueStepLength(G4double geomStepLength) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4PhysicsTable * fSecondMoments
std::vector< G4double > xsecn
G4double ComputeGeomPathLength(G4double truePathLength) override
G4ParticleChangeForMSC * fParticleChange
const G4ParticleDefinition * particle
const G4DataVector * currentCuts
const G4Material * currentMaterial
void SetupParticle(const G4ParticleDefinition *)
void StartTracking(G4Track *) override
void SetSingleScatteringFactor(G4double)
G4WentzelOKandVIxSection * wokvi
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62