Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
78G4EmDataHandler* G4GammaGeneralProcess::theHandler = nullptr;
79G4bool G4GammaGeneralProcess::theT[nTables] =
80 {true,false,true,true,true,false,true,true,true,
81 true,true,true,true,true,true};
82G4String G4GammaGeneralProcess::nameT[nTables] =
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 if (nullptr == theHandler) {
96 theHandler = new G4EmDataHandler(nTables);
97 }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
115{
116 if(nullptr == ptr) { return; }
117 G4int stype = ptr->GetProcessSubType();
118 if(stype == fRayleigh) { theRayleigh = ptr; }
119 else if(stype == fPhotoElectricEffect) { thePhotoElectric = ptr; }
120 else if(stype == fComptonScattering) { theCompton = ptr; }
121 else if(stype == fGammaConversion) { theConversionEE = ptr; }
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125
127{
128 theConversionMM = ptr;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139
141{
142 SetParticle(&part);
143 preStepLambda = 0.0;
144 idxEnergy = 0;
145 currentCouple = nullptr;
146
148
149 // initialise base material for the current run
153
154 if(1 < verboseLevel) {
155 G4cout << "G4GammaGeneralProcess::PreparePhysicsTable() for "
156 << GetProcessName()
157 << " and particle " << part.GetParticleName()
158 << " isMaster: " << isTheMaster << G4endl;
159 }
160
161 // 3 sub-processes must be always defined
162 if(thePhotoElectric == nullptr || theCompton == nullptr ||
163 theConversionEE == nullptr) {
165 ed << "### G4GeneralGammaProcess is initialized incorrectly"
166 << "\n Photoelectric: " << thePhotoElectric
167 << "\n Compton: " << theCompton
168 << "\n Conversion: " << theConversionEE;
169 G4Exception("G4GeneralGammaProcess","em0004",
170 FatalException, ed,"");
171 return;
172 }
173
174 thePhotoElectric->PreparePhysicsTable(part);
175 theCompton->PreparePhysicsTable(part);
176 theConversionEE->PreparePhysicsTable(part);
177 if (nullptr != theRayleigh) { theRayleigh->PreparePhysicsTable(part); }
178 if (nullptr != theGammaNuclear) { theGammaNuclear->PreparePhysicsTable(part); }
179 if (nullptr != theConversionMM) { theConversionMM->PreparePhysicsTable(part); }
180
181 InitialiseProcess(&part);
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187{
188 if(isTheMaster) {
189
192
193 // tables are created and its size is defined only once
194 if (nullptr != theRayleigh) { theT[1] = true; }
195
196 theHandler->SetMasterProcess(thePhotoElectric);
197 theHandler->SetMasterProcess(theCompton);
198 theHandler->SetMasterProcess(theConversionEE);
199 theHandler->SetMasterProcess(theRayleigh);
200
201 auto bld = man->GetTableBuilder();
202
203 const G4ProductionCutsTable* theCoupleTable=
205 std::size_t numOfCouples = theCoupleTable->GetTableSize();
206
207 G4double mine = param->MinKinEnergy();
208 G4double maxe = param->MaxKinEnergy();
209 G4int nd = param->NumberOfBinsPerDecade();
210 std::size_t nbin1 = std::max(5, nd*G4lrint(std::log10(minPEEnergy/mine)));
211 std::size_t nbin2 = std::max(5, nd*G4lrint(std::log10(maxe/minMMEnergy)));
212
213 G4PhysicsVector* vec = nullptr;
214 G4PhysicsLogVector aVector(mine,minPEEnergy,nbin1,true);
215 G4PhysicsLogVector bVector(minPEEnergy,minEEEnergy,nLowE,false);
216 G4PhysicsLogVector cVector(minEEEnergy,minMMEnergy,nHighE,false);
217 G4PhysicsLogVector dVector(minMMEnergy,maxe,nbin2,true);
218
219 for(std::size_t i=0; i<nTables; ++i) {
220 if(!theT[i]) { continue; }
221 //G4cout << "## PreparePhysTable " << i << "." << G4endl;
222 G4PhysicsTable* table = theHandler->MakeTable(i);
223 //G4cout << " make table " << table << G4endl;
224 for(std::size_t j=0; j<numOfCouples; ++j) {
225 vec = (*table)[j];
226 if (bld->GetFlag(j) && nullptr == vec) {
227 if(i<=1) {
228 vec = new G4PhysicsVector(aVector);
229 } else if(i<=5) {
230 vec = new G4PhysicsVector(bVector);
231 } else if(i<=9) {
232 vec = new G4PhysicsVector(cVector);
233 } else {
234 vec = new G4PhysicsVector(dVector);
235 }
237 }
238 }
239 }
240 }
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244
246{
247 if(1 < verboseLevel) {
248 G4cout << "### G4VEmProcess::BuildPhysicsTable() for "
249 << GetProcessName()
250 << " and particle " << part.GetParticleName()
251 << G4endl;
252 }
253 if(!isTheMaster) {
254 thePhotoElectric->SetEmMasterProcess(theHandler->GetMasterProcess(0));
255 baseMat = theHandler->GetMasterProcess(0)->UseBaseMaterial();
256 }
257 thePhotoElectric->BuildPhysicsTable(part);
258
259 if(!isTheMaster) {
260 theCompton->SetEmMasterProcess(theHandler->GetMasterProcess(1));
261 }
262 theCompton->BuildPhysicsTable(part);
263
264 if(!isTheMaster) {
265 theConversionEE->SetEmMasterProcess(theHandler->GetMasterProcess(2));
266 }
267 theConversionEE->BuildPhysicsTable(part);
268
269 if(theRayleigh != nullptr) {
270 if(!isTheMaster) {
271 theRayleigh->SetEmMasterProcess(theHandler->GetMasterProcess(3));
272 }
273 theRayleigh->BuildPhysicsTable(part);
274 }
275 if(theGammaNuclear != nullptr) { theGammaNuclear->BuildPhysicsTable(part); }
276 if(theConversionMM != nullptr) { theConversionMM->BuildPhysicsTable(part); }
277
278 if(isTheMaster) {
279 const G4ProductionCutsTable* theCoupleTable=
281 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
282
284 const std::vector<G4PhysicsTable*>& tables = theHandler->GetTables();
285
287 ? theGammaNuclear->GetCrossSectionDataStore() : nullptr;
288 G4DynamicParticle* dynParticle =
290
291 G4double sigComp(0.), sigPE(0.), sigConv(0.), sigR(0.),
292 sigN(0.), sigM(0.), val(0.);
293
294 for(G4int i=0; i<numOfCouples; ++i) {
295
296 if (bld->GetFlag(i)) {
297 G4int idx = (!baseMat) ? i : DensityIndex(i);
298 const G4MaterialCutsCouple* couple =
299 theCoupleTable->GetMaterialCutsCouple(i);
300 const G4Material* material = couple->GetMaterial();
301
302 // energy interval 0
303 std::size_t nn = (*(tables[0]))[idx]->GetVectorLength();
304 if(1 < verboseLevel) {
305 G4cout << "======= Zone 0 ======= N= " << nn
306 << " for " << material->GetName() << G4endl;
307 }
308 for(std::size_t j=0; j<nn; ++j) {
309 G4double e = (*(tables[0]))[idx]->Energy(j);
310 G4double loge = G4Log(e);
311 sigComp = theCompton->GetLambda(e, couple, loge);
312 sigR = (nullptr != theRayleigh) ?
313 theRayleigh->GetLambda(e, couple, loge) : 0.0;
314 G4double sum = sigComp + sigR;
315 if(1 < verboseLevel) {
316 G4cout << j << ". E= " << e << " xs= " << sum
317 << " compt= " << sigComp << " Rayl= " << sigR << G4endl;
318 }
319 (*(tables[0]))[idx]->PutValue(j, sum);
320 if(theT[1]) {
321 val = sigR/sum;
322 (*(tables[1]))[idx]->PutValue(j, val);
323 }
324 }
325
326 // energy interval 1
327 nn = (*(tables[2]))[idx]->GetVectorLength();
328 if(1 < verboseLevel) {
329 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
330 }
331 for(std::size_t j=0; j<nn; ++j) {
332 G4double e = (*(tables[2]))[idx]->Energy(j);
333 G4double loge = G4Log(e);
334 sigComp = theCompton->GetLambda(e, couple, loge);
335 sigR = (nullptr != theRayleigh) ?
336 theRayleigh->GetLambda(e, couple, loge) : 0.0;
337 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
338 G4double sum = sigComp + sigR + sigPE;
339 if(1 < verboseLevel) {
340 G4cout << j << ". E= " << e << " xs= " << sum
341 << " compt= " << sigComp << " conv= " << sigConv
342 << " PE= " << sigPE << " Rayl= " << sigR
343 << " GN= " << sigN << G4endl;
344 }
345 (*(tables[2]))[idx]->PutValue(j, sum);
346
347 val = sigPE/sum;
348 (*(tables[3]))[idx]->PutValue(j, val);
349
350 val = (sigR > 0.0) ? (sigComp + sigPE)/sum : 1.0;
351 (*(tables[4]))[idx]->PutValue(j, val);
352 }
353
354 // energy interval 2
355 nn = (*(tables[6]))[idx]->GetVectorLength();
356 if(1 < verboseLevel) {
357 G4cout << "======= Zone 2 ======= N= " << nn << G4endl;
358 }
359 for(std::size_t j=0; j<nn; ++j) {
360 G4double e = (*(tables[6]))[idx]->Energy(j);
361 G4double loge = G4Log(e);
362 sigComp = theCompton->GetLambda(e, couple, loge);
363 sigConv = theConversionEE->GetLambda(e, couple, loge);
364 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
365 sigN = 0.0;
366 if(nullptr != gn) {
367 dynParticle->SetKineticEnergy(e);
368 sigN = gn->ComputeCrossSection(dynParticle, material);
369 }
370 G4double sum = sigComp + sigConv + sigPE + sigN;
371 if(1 < verboseLevel) {
372 G4cout << j << ". E= " << e << " xs= " << sum
373 << " compt= " << sigComp << " conv= " << sigConv
374 << " PE= " << sigPE
375 << " GN= " << sigN << G4endl;
376 }
377 (*(tables[6]))[idx]->PutValue(j, sum);
378
379 val = sigConv/sum;
380 (*(tables[7]))[idx]->PutValue(j, val);
381
382 val = (sigConv + sigComp)/sum;
383 (*(tables[8]))[idx]->PutValue(j, val);
384
385 val = (sigN > 0.0) ? (sigConv + sigComp + sigPE)/sum : 1.0;
386 (*(tables[9]))[idx]->PutValue(j, val);
387 }
388
389 // energy interval 3
390 nn = (*(tables[10]))[idx]->GetVectorLength();
391 if(1 < verboseLevel) {
392 G4cout << "======= Zone 3 ======= N= " << nn
393 << " for " << material->GetName() << G4endl;
394 }
395 for(std::size_t j=0; j<nn; ++j) {
396 G4double e = (*(tables[10]))[idx]->Energy(j);
397 G4double loge = G4Log(e);
398 sigComp = theCompton->GetLambda(e, couple, loge);
399 sigConv = theConversionEE->GetLambda(e, couple, loge);
400 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
401 sigN = 0.0;
402 if(nullptr != gn) {
403 dynParticle->SetKineticEnergy(e);
404 sigN = gn->ComputeCrossSection(dynParticle, material);
405 }
406 sigM = 0.0;
407 if(nullptr != theConversionMM) {
408 val = theConversionMM->ComputeMeanFreePath(e, material);
409 sigM = (val < DBL_MAX) ? 1./val : 0.0;
410 }
411 G4double sum = sigComp + sigConv + sigPE + sigN + sigM;
412 if(1 < verboseLevel) {
413 G4cout << j << ". E= " << e << " xs= " << sum
414 << " compt= " << sigComp << " conv= " << sigConv
415 << " PE= " << sigPE
416 << " GN= " << sigN << G4endl;
417 }
418 (*(tables[10]))[idx]->PutValue(j, sum);
419
420 val = (sigComp + sigPE + sigN + sigM)/sum;
421 (*(tables[11]))[idx]->PutValue(j, val);
422
423 val = (sigPE + sigN + sigM)/sum;
424 (*(tables[12]))[idx]->PutValue(j, val);
425
426 val = (sigN + sigM)/sum;
427 (*(tables[13]))[idx]->PutValue(j, val);
428
429 val = sigM/sum;
430 (*(tables[14]))[idx]->PutValue(j, val);
431 }
432 for(std::size_t k=0; k<nTables; ++k) {
433 if(theT[k] && (k <= 1 || k >= 10)) {
434 //G4cout <<"BuildPhysTable spline iTable="<<k<<" jCouple="<< idx << G4endl;
435 (*(tables[k]))[idx]->FillSecondDerivatives();
436 }
437 }
438 }
439 }
440 delete dynParticle;
441 }
442
443 if(1 < verboseLevel) {
444 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
445 << GetProcessName()
446 << " and particle " << part.GetParticleName()
447 << G4endl;
448 }
449}
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
452
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459
461 const G4Track& track,
462 G4double previousStepSize,
464{
466 G4double x = DBL_MAX;
467
468 G4double energy = track.GetKineticEnergy();
469 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
470
471 // compute mean free path
472 G4bool recompute = false;
473 if(couple != currentCouple) {
474 currentCouple = couple;
476 currentMaterial = couple->GetMaterial();
477 factor = 1.0;
478 if(baseMat) {
481 }
482 recompute = true;
483 }
484 if(energy != preStepKinEnergy) {
485 preStepKinEnergy = energy;
487 recompute = true;
488 }
489 if(recompute) {
491
492 // zero cross section
493 if(preStepLambda <= 0.0) {
496 }
497 }
498
499 // non-zero cross section
500 if(preStepLambda > 0.0) {
501
503
504 // beggining of tracking (or just after DoIt of this process)
507
508 } else if(currentInteractionLength < DBL_MAX) {
509
511 previousStepSize/currentInteractionLength;
514 }
515
516 // new mean free path and step limit for the next step
519 }
520 /*
521 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
522 << " idxe= " << idxEnergy << " xs= " << preStepLambda
523 << " x= " << x << G4endl;
524 */
525 return x;
526}
527
528//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
529
531{
532 G4double cross = 0.0;
533 /*
534 G4cout << "#Total: " << preStepKinEnergy << " " << minPEEnergy << " "
535 << minEEEnergy << " " << minMMEnergy<< G4endl;
536 G4cout << " idxE= " << idxEnergy
537 << " idxC= " << currentCoupleIndex << G4endl;
538 */
539 if(preStepKinEnergy < minPEEnergy) {
540 cross = ComputeGeneralLambda(0, 0);
541 //G4cout << "XS1: " << cross << G4endl;
542 peLambda = thePhotoElectric->GetLambda(preStepKinEnergy, currentCouple, preStepLogE);
543 cross += peLambda;
544 //G4cout << "XS2: " << peLambda << G4endl;
545
546 } else if(preStepKinEnergy < minEEEnergy) {
547 cross = ComputeGeneralLambda(1, 2);
548 //G4cout << "XS3: " << cross << G4endl;
549
550 } else if(preStepKinEnergy < minMMEnergy) {
551 cross = ComputeGeneralLambda(2, 6);
552 //G4cout << "XS4: " << cross << G4endl;
553
554 } else {
555 cross = ComputeGeneralLambda(3, 10);
556 //G4cout << "XS5: " << cross << G4endl;
557 }
558 /*
559 G4cout << "xs= " << cross << " idxE= " << idxEnergy
560 << " idxC= " << currentCoupleIndex
561 << " E= " << preStepKinEnergy << G4endl;
562 */
563 return cross;
564}
565
566//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
567
569 const G4Step& step)
570{
571 // In all cases clear number of interaction lengths
573 selectedProc = nullptr;
575 /*
576 G4cout << "PostStep: preStepLambda= " << preStepLambda
577 << " PE= " << peLambda << " q= " << q << " idxE= " << idxEnergy
578 << G4endl;
579 */
580 switch (idxEnergy) {
581 case 0:
582 if(preStepLambda*q <= peLambda) {
583 SelectEmProcess(step, thePhotoElectric);
584 } else {
585 if(theT[1] && preStepLambda*q < (preStepLambda - peLambda)*GetProbability(1) + peLambda) {
586 SelectEmProcess(step, theRayleigh);
587 } else {
588 SelectEmProcess(step, theCompton);
589 }
590 }
591 break;
592
593 case 1:
594 if(q <= GetProbability(3)) {
595 SelectEmProcess(step, thePhotoElectric);
596 } else if(q <= GetProbability(4)) {
597 SelectEmProcess(step, theCompton);
598 } else if(theRayleigh) {
599 SelectEmProcess(step, theRayleigh);
600 } else {
601 SelectEmProcess(step, thePhotoElectric);
602 }
603 break;
604
605 case 2:
606 if(q <= GetProbability(7)) {
607 SelectEmProcess(step, theConversionEE);
608 } else if(q <= GetProbability(8)) {
609 SelectEmProcess(step, theCompton);
610 } else if(q <= GetProbability(9)) {
611 SelectEmProcess(step, thePhotoElectric);
612 } else if(theGammaNuclear) {
613 SelectHadProcess(track, step, theGammaNuclear);
614 } else {
615 SelectEmProcess(step, theConversionEE);
616 }
617 break;
618
619 case 3:
620 if(q + GetProbability(11) <= 1.0) {
621 SelectEmProcess(step, theConversionEE);
622 } else if(q + GetProbability(12) <= 1.0) {
623 SelectEmProcess(step, theCompton);
624 } else if(q + GetProbability(13) <= 1.0) {
625 SelectEmProcess(step, thePhotoElectric);
626 } else if(theGammaNuclear && q + GetProbability(14) <= 1.0) {
627 SelectHadProcess(track, step, theGammaNuclear);
628 } else if(theConversionMM) {
629 SelectedProcess(step, theConversionMM);
630 } else {
631 SelectEmProcess(step, theConversionEE);
632 }
633 break;
634 }
635 // sample secondaries
636 if(selectedProc != nullptr) {
637 return selectedProc->PostStepDoIt(track, step);
638 }
639 // no interaction - exception case
640 fParticleChange.InitializeForPostStep(track);
641 return &fParticleChange;
642}
643
644//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
645
653
654//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
655
657 const G4String& directory,
658 G4bool ascii)
659{
660 G4bool yes = true;
661 if(!isTheMaster) { return yes; }
662 if(!thePhotoElectric->StorePhysicsTable(part, directory, ascii))
663 { yes = false; }
664 if(!theCompton->StorePhysicsTable(part, directory, ascii))
665 { yes = false; }
666 if(!theConversionEE->StorePhysicsTable(part, directory, ascii))
667 { yes = false; }
668 if(theRayleigh != nullptr &&
669 !theRayleigh->StorePhysicsTable(part, directory, ascii))
670 { yes = false; }
671
672 for(std::size_t i=0; i<nTables; ++i) {
673 if(theT[i]) {
674 G4String nam = (0==i || 2==i || 6==i || 10==i)
675 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
676 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
677 if(!theHandler->StorePhysicsTable(i, part, fnam, ascii)) { yes = false; }
678 }
679 }
680 return yes;
681}
682
683//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
684
685G4bool
687 const G4String& directory,
688 G4bool ascii)
689{
690 if (!isTheMaster) { return true; }
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(std::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 G4bool spline = (i <= 1 || i >= 10);
713 if(!theHandler->RetrievePhysicsTable(i, part, fnam, ascii, spline))
714 { yes = false; }
715 }
716 }
717 return yes;
718}
719
720//....Ooooo0ooooo ........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
721
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
732void G4GammaGeneralProcess::ProcessDescription(std::ostream& out) const
733{
734 thePhotoElectric->ProcessDescription(out);
735 theCompton->ProcessDescription(out);
736 theConversionEE->ProcessDescription(out);
737 if(theRayleigh) { theRayleigh->ProcessDescription(out); }
738 if(theGammaNuclear) { theGammaNuclear->ProcessDescription(out); }
739 if(theConversionMM) { theConversionMM->ProcessDescription(out); }
740}
741
742//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
743
745{
746 return (selectedProc) ? selectedProc->GetProcessName()
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
751
753{
754 return (selectedProc) ? selectedProc->GetProcessSubType()
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
761{
762 G4VEmProcess* proc = nullptr;
763 if(name == thePhotoElectric->GetProcessName()) {
764 proc = thePhotoElectric;
765 } else if(name == theCompton->GetProcessName()) {
766 proc = theCompton;
767 } else if(name == theConversionEE->GetProcessName()) {
768 proc = theConversionEE;
769 } else if(theRayleigh != nullptr && name == theRayleigh->GetProcessName()) {
770 proc = theRayleigh;
771 }
772 return proc;
773}
774
775//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
776
781
782//....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)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition G4Log.hh:169
@ fElectromagnetic
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:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetLogKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
void AddHadProcess(G4HadronicProcess *)
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
G4bool IsApplicable(const G4ParticleDefinition &) override
void AddEmProcess(G4VEmProcess *)
void StartTracking(G4Track *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void AddMMProcess(G4GammaConversionToMuons *)
void ProcessDescription(std::ostream &outFile) const override
G4double ComputeGeneralLambda(std::size_t idxe, std::size_t idxt)
const G4String & GetSubProcessName() const
G4GammaGeneralProcess(const G4String &pname="GammaGeneralProc")
G4VEmProcess * GetEmProcess(const G4String &name) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
const G4VProcess * GetCreatorProcess() const override
void InitialiseProcess(const G4ParticleDefinition *) override
void SelectEmProcess(const G4Step &, G4VEmProcess *)
G4double GetProbability(std::size_t idxt)
G4HadronicProcess * theGammaNuclear
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:81
G4CrossSectionDataStore * GetCrossSectionDataStore()
static G4bool GetBaseMaterialFlag()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
static G4bool GetFlag(std::size_t idx)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
const G4Material * GetMaterial() const
const G4String & GetName() const
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()
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double MeanFreePath(const G4Track &track)
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
G4double preStepLambda
G4int DensityIndex(G4int idx) const
std::size_t currentCoupleIndex
std::size_t basedCoupleIndex
G4double DensityFactor(G4int idx) const
const G4MaterialCutsCouple * currentCouple
G4double preStepKinEnergy
G4ParticleChangeForGamma fParticleChange
void SetParticle(const G4ParticleDefinition *p)
const G4Material * currentMaterial
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength
void SetVerboseLevel(G4int value)
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46
G4int verboseLevel
G4double theNumberOfInteractionLengthLeft
void SetProcessSubType(G4int)
G4int GetProcessSubType() const
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
const G4String & GetProcessName() const
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62