Geant4-11
G4GammaGeneralProcess.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: G4GammaGeneralProcess
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 19.07.2018
37//
38// Modifications:
39//
40// Class Description:
41//
42
43// -------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
49#include "G4VDiscreteProcess.hh"
51#include "G4SystemOfUnits.hh"
52#include "G4ProcessManager.hh"
54#include "G4LossTableBuilder.hh"
55#include "G4HadronicProcess.hh"
56#include "G4LossTableManager.hh"
57#include "G4Step.hh"
58#include "G4Track.hh"
60#include "G4DataVector.hh"
61#include "G4PhysicsTable.hh"
62#include "G4PhysicsLogVector.hh"
64#include "G4VParticleChange.hh"
66#include "G4EmParameters.hh"
67#include "G4EmProcessSubType.hh"
68#include "G4Material.hh"
71#include "G4Gamma.hh"
72
73#include "G4Log.hh"
74#include <iostream>
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
80 {true,false,true,true,true,false,true,true,true,
81 true,true,true,true,true,true};
83 {"0","1","2","3","4","5","6","7","8",
84 "9","10","11","12","13","14"};
85
88 minPEEnergy(150*CLHEP::keV),
89 minEEEnergy(2*CLHEP::electron_mass_c2),
90 minMMEnergy(100*CLHEP::MeV)
91{
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
100{
101 if(isTheMaster) {
102 delete theHandler;
103 theHandler = nullptr;
104 }
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108
110{
111 return true;
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
117{
118 if(nullptr == ptr) { return; }
119 G4int stype = ptr->GetProcessSubType();
120 if(stype == fRayleigh) { theRayleigh = ptr; }
121 else if(stype == fPhotoElectricEffect) { thePhotoElectric = ptr; }
122 else if(stype == fComptonScattering) { theCompton = ptr; }
123 else if(stype == fGammaConversion) { theConversionEE = ptr; }
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
129{
130 theConversionMM = ptr;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134
136{
137 theGammaNuclear = ptr;
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
143{
144 SetParticle(&part);
145 preStepLambda = 0.0;
146 idxEnergy = 0;
147 currentCouple = nullptr;
148
151
152 isTheMaster = man->IsMaster();
153 if(isTheMaster) { SetVerboseLevel(param->Verbose()); }
154 else { SetVerboseLevel(param->WorkerVerbose()); }
155
158
159 if(1 < verboseLevel) {
160 G4cout << "G4GammaGeneralProcess::PreparePhysicsTable() for "
161 << GetProcessName()
162 << " and particle " << part.GetParticleName()
163 << " isMaster: " << isTheMaster << G4endl;
164 }
165
166 // 3 sub-processes must be always defined
167 if(thePhotoElectric == nullptr || theCompton == nullptr ||
168 theConversionEE == nullptr) {
170 ed << "### G4GeneralGammaProcess is initialized incorrectly"
171 << "\n Photoelectric: " << thePhotoElectric
172 << "\n Compton: " << theCompton
173 << "\n Conversion: " << theConversionEE;
174 G4Exception("G4GeneralGammaProcess","em0004",
175 FatalException, ed,"");
176 }
177
184
185 InitialiseProcess(&part);
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
191{
192 if(isTheMaster) {
193
196
197 // tables are created and its size is defined only once
198 if(nullptr == theHandler) {
200 if(theRayleigh) { theT[1] = true; }
201
206 }
207 auto bld = man->GetTableBuilder();
208
209 const G4ProductionCutsTable* theCoupleTable=
211 size_t numOfCouples = theCoupleTable->GetTableSize();
212
213 G4double mine = param->MinKinEnergy();
214 G4double maxe = param->MaxKinEnergy();
215 G4int nd = param->NumberOfBinsPerDecade();
216 size_t nbin1 = std::max(5, nd*G4lrint(std::log10(minPEEnergy/mine)));
217 size_t nbin2 = std::max(5, nd*G4lrint(std::log10(maxe/minMMEnergy)));
218
219 G4PhysicsVector* vec = nullptr;
220 G4PhysicsLogVector aVector(mine,minPEEnergy,nbin1,splineFlag);
223 G4PhysicsLogVector dVector(minMMEnergy,maxe,nbin2,splineFlag);
224
225 for(size_t i=0; i<nTables; ++i) {
226 if(!theT[i]) { continue; }
227 //G4cout << "## PreparePhysTable " << i << "." << G4endl;
229 //G4cout << " make table " << table << G4endl;
230 for(size_t j=0; j<numOfCouples; ++j) {
231 vec = (*table)[j];
232 if (bld->GetFlag(j) && nullptr == vec) {
233 //G4cout <<" i= "<<i<<" j= "<< j <<" make new vector"<< G4endl;
234 if(i<=1) {
235 vec = new G4PhysicsVector(aVector);
236 } else if(i<=5) {
237 vec = new G4PhysicsVector(bVector);
238 } else if(i<=9) {
239 vec = new G4PhysicsVector(cVector);
240 } else {
241 vec = new G4PhysicsVector(dVector);
242 }
244 }
245 }
246 }
247 }
248}
249
250//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
251
253{
254 if(1 < verboseLevel) {
255 G4cout << "### G4VEmProcess::BuildPhysicsTable() for "
256 << GetProcessName()
257 << " and particle " << part.GetParticleName()
258 << G4endl;
259 }
260 if(!isTheMaster) {
263 }
265
266 if(!isTheMaster) {
268 }
270
271 if(!isTheMaster) {
273 }
275
276 if(theRayleigh != nullptr) {
277 if(!isTheMaster) {
279 }
281 }
282 if(theGammaNuclear != nullptr) { theGammaNuclear->BuildPhysicsTable(part); }
283 if(theConversionMM != nullptr) { theConversionMM->BuildPhysicsTable(part); }
284
285 if(isTheMaster) {
286 const G4ProductionCutsTable* theCoupleTable=
288 size_t numOfCouples = theCoupleTable->GetTableSize();
289
291 const std::vector<G4PhysicsTable*>& tables = theHandler->GetTables();
292
295 G4DynamicParticle* dynParticle =
297
298 G4double sigComp(0.), sigPE(0.), sigConv(0.), sigR(0.),
299 sigN(0.), sigM(0.), val(0.);
300
301 for(size_t i=0; i<numOfCouples; ++i) {
302
303 if (bld->GetFlag(i)) {
304
305 G4int idx = (!baseMat) ? i : DensityIndex(i);
306 const G4MaterialCutsCouple* couple =
307 theCoupleTable->GetMaterialCutsCouple(i);
308 const G4Material* material = couple->GetMaterial();
309
310 // energy interval 0
311 size_t nn = (*(tables[0]))[idx]->GetVectorLength();
312 if(1 < verboseLevel) {
313 G4cout << "======= Zone 0 ======= N= " << nn
314 << " for " << material->GetName() << G4endl;
315 }
316 for(size_t j=0; j<nn; ++j) {
317 G4double e = (*(tables[0]))[idx]->Energy(j);
318 G4double loge = G4Log(e);
319 sigComp = theCompton->GetLambda(e, couple, loge);
320 sigR = (nullptr != theRayleigh) ?
321 theRayleigh->GetLambda(e, couple, loge) : 0.0;
322 G4double sum = sigComp + sigR;
323 if(1 < verboseLevel) {
324 G4cout << j << ". E= " << e << " xs= " << sum
325 << " compt= " << sigComp << " Rayl= " << sigR << G4endl;
326 }
327 (*(tables[0]))[idx]->PutValue(j, sum);
328 if(theT[1]) {
329 val = sigR/sum;
330 (*(tables[1]))[idx]->PutValue(j, val);
331 }
332 }
333
334 // energy interval 1
335 nn = (*(tables[2]))[idx]->GetVectorLength();
336 if(1 < verboseLevel) {
337 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
338 }
339 for(size_t j=0; j<nn; ++j) {
340 G4double e = (*(tables[2]))[idx]->Energy(j);
341 G4double loge = G4Log(e);
342 sigComp = theCompton->GetLambda(e, couple, loge);
343 sigR = (nullptr != theRayleigh) ?
344 theRayleigh->GetLambda(e, couple, loge) : 0.0;
345 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
346 G4double sum = sigComp + sigR + sigPE;
347 if(1 < verboseLevel) {
348 G4cout << j << ". E= " << e << " xs= " << sum
349 << " compt= " << sigComp << " conv= " << sigConv
350 << " PE= " << sigPE << " Rayl= " << sigR
351 << " GN= " << sigN << G4endl;
352 }
353 (*(tables[2]))[idx]->PutValue(j, sum);
354
355 val = sigPE/sum;
356 (*(tables[3]))[idx]->PutValue(j, val);
357
358 val = (sigR > 0.0) ? (sigComp + sigPE)/sum : 1.0;
359 (*(tables[4]))[idx]->PutValue(j, val);
360 }
361
362 // energy interval 2
363 nn = (*(tables[6]))[idx]->GetVectorLength();
364 if(1 < verboseLevel) {
365 G4cout << "======= Zone 2 ======= N= " << nn << G4endl;
366 }
367 for(size_t j=0; j<nn; ++j) {
368 G4double e = (*(tables[6]))[idx]->Energy(j);
369 G4double loge = G4Log(e);
370 sigComp = theCompton->GetLambda(e, couple, loge);
371 sigConv = theConversionEE->GetLambda(e, couple, loge);
372 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
373 sigN = 0.0;
374 if(nullptr != gn) {
375 dynParticle->SetKineticEnergy(e);
376 sigN = gn->ComputeCrossSection(dynParticle, material);
377 }
378 G4double sum = sigComp + sigConv + sigPE + sigN;
379 if(1 < verboseLevel) {
380 G4cout << j << ". E= " << e << " xs= " << sum
381 << " compt= " << sigComp << " conv= " << sigConv
382 << " PE= " << sigPE
383 << " GN= " << sigN << G4endl;
384 }
385 (*(tables[6]))[idx]->PutValue(j, sum);
386
387 val = sigConv/sum;
388 (*(tables[7]))[idx]->PutValue(j, val);
389
390 val = (sigConv + sigComp)/sum;
391 (*(tables[8]))[idx]->PutValue(j, val);
392
393 val = (sigN > 0.0) ? (sigConv + sigComp + sigPE)/sum : 1.0;
394 (*(tables[9]))[idx]->PutValue(j, val);
395 }
396
397 // energy interval 3
398 nn = (*(tables[10]))[idx]->GetVectorLength();
399 if(1 < verboseLevel) {
400 G4cout << "======= Zone 3 ======= N= " << nn
401 << " for " << material->GetName() << G4endl;
402 }
403 for(size_t j=0; j<nn; ++j) {
404 G4double e = (*(tables[10]))[idx]->Energy(j);
405 G4double loge = G4Log(e);
406 sigComp = theCompton->GetLambda(e, couple, loge);
407 sigConv = theConversionEE->GetLambda(e, couple, loge);
408 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
409 sigN = 0.0;
410 if(nullptr != gn) {
411 dynParticle->SetKineticEnergy(e);
412 sigN = gn->ComputeCrossSection(dynParticle, material);
413 }
414 sigM = 0.0;
415 if(nullptr != theConversionMM) {
417 sigM = (val < DBL_MAX) ? 1./val : 0.0;
418 }
419 G4double sum = sigComp + sigConv + sigPE + sigN + sigM;
420 if(1 < verboseLevel) {
421 G4cout << j << ". E= " << e << " xs= " << sum
422 << " compt= " << sigComp << " conv= " << sigConv
423 << " PE= " << sigPE
424 << " GN= " << sigN << G4endl;
425 }
426 (*(tables[10]))[idx]->PutValue(j, sum);
427
428 val = (sigComp + sigPE + sigN + sigM)/sum;
429 (*(tables[11]))[idx]->PutValue(j, val);
430
431 val = (sigPE + sigN + sigM)/sum;
432 (*(tables[12]))[idx]->PutValue(j, val);
433
434 val = (sigN + sigM)/sum;
435 (*(tables[13]))[idx]->PutValue(j, val);
436
437 val = sigN/sum;
438 (*(tables[14]))[idx]->PutValue(j, val);
439 }
440 for(size_t k=0; k<nTables; ++k) {
441 if(splineFlag) {
442 (*(tables[k]))[idx]->FillSecondDerivatives();
443 }
444 }
445 }
446 }
447 delete dynParticle;
448 }
449
450 if(1 < verboseLevel) {
451 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
452 << GetProcessName()
453 << " and particle " << part.GetParticleName()
454 << G4endl;
455 }
456}
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459
461{
463}
464
465//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466
468 const G4Track& track,
469 G4double previousStepSize,
471{
473 G4double x = DBL_MAX;
474
476 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
477
478 // compute mean free path
479 G4bool recompute = false;
480 if(couple != currentCouple) {
481 currentCouple = couple;
483 currentMaterial = couple->GetMaterial();
484 factor = 1.0;
485 if(baseMat) {
488 }
489 recompute = true;
490 }
491 if(energy != preStepKinEnergy) {
494 recompute = true;
495 }
496 if(recompute) {
498
499 // zero cross section
500 if(preStepLambda <= 0.0) {
503 }
504 }
505
506 // non-zero cross section
507 if(preStepLambda > 0.0) {
508
510
511 // beggining of tracking (or just after DoIt of this process)
514
515 } else if(currentInteractionLength < DBL_MAX) {
516
518 previousStepSize/currentInteractionLength;
521 }
522
523 // new mean free path and step limit for the next step
526 }
527 /*
528 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
529 << " idxe= " << idxEnergy << " xs= " << preStepLambda
530 << " x= " << x << G4endl;
531 */
532 return x;
533}
534
535//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
536
538{
539 G4double cross = 0.0;
540 /*
541 G4cout << "#Total: " << preStepKinEnergy << " " << minPEEnergy << " "
542 << minEEEnergy << " " << minMMEnergy<< G4endl;
543 G4cout << " idxE= " << idxEnergy
544 << " idxC= " << currentCoupleIndex << G4endl;
545 */
547 cross = ComputeGeneralLambda(0, 0);
548 //G4cout << "XS1: " << cross << G4endl;
550 cross += peLambda;
551 //G4cout << "XS2: " << cross << G4endl;
552
553 } else if(preStepKinEnergy < minEEEnergy) {
554 cross = ComputeGeneralLambda(1, 2);
555 //G4cout << "XS3: " << cross << G4endl;
556
557 } else if(preStepKinEnergy < minMMEnergy) {
558 cross = ComputeGeneralLambda(2, 6);
559 //G4cout << "XS4: " << cross << G4endl;
560
561 } else {
562 cross = ComputeGeneralLambda(3, 10);
563 //G4cout << "XS5: " << cross << G4endl;
564 }
565 /*
566 G4cout << "xs= " << cross << " idxE= " << idxEnergy
567 << " idxC= " << currentCoupleIndex
568 << " E= " << energy << G4endl;
569 */
570 return cross;
571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
574
576 const G4Step& step)
577{
578 // In all cases clear number of interaction lengths
580 selectedProc = nullptr;
582 /*
583 G4cout << "PostStep: preStepLambda= " << preStepLambda
584 << " PE= " << peLambda << " q= " << q << " idxE= " << idxEnergy
585 << G4endl;
586 */
587 switch (idxEnergy) {
588 case 0:
589 if(preStepLambda*q <= peLambda) {
591 } else {
594 } else {
596 }
597 }
598 break;
599
600 case 1:
601 if(q <= GetProbability(3)) {
603 } else if(q <= GetProbability(4)) {
605 } else if(theRayleigh) {
607 }
608 break;
609
610 case 2:
611 if(q <= GetProbability(7)) {
613 } else if(q <= GetProbability(8)) {
615 } else if(q <= GetProbability(9)) {
617 } else if(theGammaNuclear) {
618 SelectHadProcess(track, step, theGammaNuclear);
619 }
620 break;
621
622 case 3:
623 if(q + GetProbability(11) <= 1.0) {
625 } else if(q + GetProbability(12) <= 1.0) {
627 } else if(q + GetProbability(13) <= 1.0) {
629 } else if(theGammaNuclear && q + GetProbability(14) <= 1.0) {
630 SelectHadProcess(track, step, theGammaNuclear);
631 } else if(theConversionMM) {
633 }
634 break;
635 }
636 // sample secondaries
637 if(selectedProc != nullptr) {
638 return selectedProc->PostStepDoIt(track, step);
639 }
640 // no interaction - exception case
642 return &fParticleChange;
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646
648 const G4Step& step, G4HadronicProcess* proc)
649{
650 SelectedProcess(step, proc);
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
658 const G4String& directory,
659 G4bool ascii)
660{
661 G4bool yes = true;
662 if(!isTheMaster) { return yes; }
663 if(!thePhotoElectric->StorePhysicsTable(part, directory, ascii))
664 { yes = false; }
665 if(!theCompton->StorePhysicsTable(part, directory, ascii))
666 { yes = false; }
667 if(!theConversionEE->StorePhysicsTable(part, directory, ascii))
668 { yes = false; }
669 if(theRayleigh != nullptr &&
670 !theRayleigh->StorePhysicsTable(part, directory, ascii))
671 { yes = false; }
672
673 for(size_t i=0; i<nTables; ++i) {
674 if(theT[i]) {
675 G4String nam = (0==i || 2==i || 6==i || 10==i)
676 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
677 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
678 if(!theHandler->StorePhysicsTable(i, part, fnam, ascii)) { yes = false; }
679 }
680 }
681 return yes;
682}
683
684//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
685
686G4bool
688 const G4String& directory,
689 G4bool ascii)
690{
691 if(1 < verboseLevel) {
692 G4cout << "G4GammaGeneralProcess::RetrievePhysicsTable() for "
693 << part->GetParticleName() << " and process "
694 << GetProcessName() << G4endl;
695 }
696 G4bool yes = true;
697 if(!thePhotoElectric->RetrievePhysicsTable(part, directory, ascii))
698 { yes = false; }
699 if(!theCompton->RetrievePhysicsTable(part, directory, ascii))
700 { yes = false; }
701 if(!theConversionEE->RetrievePhysicsTable(part, directory, ascii))
702 { yes = false; }
703 if(theRayleigh != nullptr &&
704 !theRayleigh->RetrievePhysicsTable(part, directory, ascii))
705 { yes = false; }
706
707 for(size_t i=0; i<nTables; ++i) {
708 if(theT[i]) {
709 G4String nam = (0==i || 2==i || 6==i || 10==i)
710 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
711 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
712 if(!theHandler->RetrievePhysicsTable(i, part, fnam, ascii, splineFlag))
713 { yes = false; }
714 }
715 }
716 if(yes) {
717 }
718 return yes;
719}
720
721//....Ooooo0ooooo ........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
722
724 G4double,
726{
728 return MeanFreePath(track);
729}
730
731//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
732
733void G4GammaGeneralProcess::ProcessDescription(std::ostream& out) const
734{
741}
742
743//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
744
746{
749}
750
751//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752
754{
757}
758
759//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
760
762{
763 G4VEmProcess* proc = nullptr;
765 proc = thePhotoElectric;
766 } else if(name == theCompton->GetProcessName()) {
767 proc = theCompton;
768 } else if(name == theConversionEE->GetProcessName()) {
769 proc = theConversionEE;
770 } else if(theRayleigh != nullptr && name == theRayleigh->GetProcessName()) {
771 proc = theRayleigh;
772 }
773 return proc;
774}
775
776//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fGammaGeneralProcess
@ fGammaConversion
@ fRayleigh
@ fComptonScattering
@ fPhotoElectricEffect
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fElectromagnetic
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double MeV
Definition: G4SIunits.hh:200
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetLogKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
void SetMasterProcess(const G4VEmProcess *)
G4PhysicsTable * MakeTable(size_t idx)
G4bool RetrievePhysicsTable(size_t idx, const G4ParticleDefinition *part, const G4String &fname, G4bool ascii, G4bool spline)
const G4VEmProcess * GetMasterProcess(size_t idx) const
const std::vector< G4PhysicsTable * > & GetTables() const
G4bool StorePhysicsTable(size_t idx, const G4ParticleDefinition *part, const G4String &fname, G4bool ascii)
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4int Verbose() const
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
G4double ComputeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
void BuildPhysicsTable(const G4ParticleDefinition &) override
void AddHadProcess(G4HadronicProcess *)
static G4String nameT[nTables]
void SelectedProcess(const G4Step &step, G4VProcess *ptr)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void SelectHadProcess(const G4Track &, const G4Step &, G4HadronicProcess *)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetProbability(size_t idxt)
G4bool IsApplicable(const G4ParticleDefinition &) override
void AddEmProcess(G4VEmProcess *)
void StartTracking(G4Track *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void AddMMProcess(G4GammaConversionToMuons *)
G4double ComputeGeneralLambda(size_t idxe, size_t idxt)
void ProcessDescription(std::ostream &outFile) const override
const G4String & GetSubProcessName() const
G4GammaGeneralProcess(const G4String &pname="GammaGeneralProc")
G4VEmProcess * GetEmProcess(const G4String &name) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
void InitialiseProcess(const G4ParticleDefinition *) override
void SelectEmProcess(const G4Step &, G4VEmProcess *)
static G4EmDataHandler * theHandler
G4HadronicProcess * theGammaNuclear
static const size_t nTables
G4GammaConversionToMuons * theConversionMM
static G4bool theT[nTables]
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void ProcessDescription(std::ostream &outFile) const override
void PreparePhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4bool GetFlag(size_t idx)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
const G4Material * GetMaterial() const
void InitializeForPostStep(const G4Track &)
const G4String & GetParticleName() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
Definition: G4Step.hh:62
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double MeanFreePath(const G4Track &track)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
size_t basedCoupleIndex
G4int DensityIndex(G4int idx) const
void SetEmMasterProcess(const G4VEmProcess *)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4bool UseBaseMaterial() const
void ProcessDescription(std::ostream &outFile) const override
G4double DensityFactor(G4int idx) const
G4bool isTheMaster
const G4MaterialCutsCouple * currentCouple
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double preStepKinEnergy
size_t currentCoupleIndex
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4ParticleChangeForGamma fParticleChange
void SetParticle(const G4ParticleDefinition *p)
void PreparePhysicsTable(const G4ParticleDefinition &) override
const G4Material * currentMaterial
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:175
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
Definition: DoubConv.h:17
G4double energy(const ThreeVector &p, const G4double m)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
string material
Definition: eplot.py:19
string pname
Definition: eplot.py:33
float electron_mass_c2
Definition: hepunit.py:273
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62