Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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}
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
178 if(thePhotoElectric) { thePhotoElectric->PreparePhysicsTable(part); }
179 if(theCompton) { theCompton->PreparePhysicsTable(part); }
180 if(theConversionEE) { theConversionEE->PreparePhysicsTable(part); }
181 if(theRayleigh) { theRayleigh->PreparePhysicsTable(part); }
183 if(theConversionMM) { theConversionMM->PreparePhysicsTable(part); }
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) {
199 theHandler = new G4EmDataHandler(nTables);
200 if(theRayleigh) { theT[1] = true; }
201
202 theHandler->SetMasterProcess(thePhotoElectric);
203 theHandler->SetMasterProcess(theCompton);
204 theHandler->SetMasterProcess(theConversionEE);
205 theHandler->SetMasterProcess(theRayleigh);
206 }
207 auto bld = man->GetTableBuilder();
208
209 const G4ProductionCutsTable* theCoupleTable=
211 std::size_t numOfCouples = theCoupleTable->GetTableSize();
212
213 G4double mine = param->MinKinEnergy();
214 G4double maxe = param->MaxKinEnergy();
215 G4int nd = param->NumberOfBinsPerDecade();
216 std::size_t nbin1 = std::max(5, nd*G4lrint(std::log10(minPEEnergy/mine)));
217 std::size_t nbin2 = std::max(5, nd*G4lrint(std::log10(maxe/minMMEnergy)));
218
219 G4PhysicsVector* vec = nullptr;
220 G4PhysicsLogVector aVector(mine,minPEEnergy,nbin1,true);
221 G4PhysicsLogVector bVector(minPEEnergy,minEEEnergy,nLowE,false);
222 G4PhysicsLogVector cVector(minEEEnergy,minMMEnergy,nHighE,false);
223 G4PhysicsLogVector dVector(minMMEnergy,maxe,nbin2,true);
224
225 for(std::size_t i=0; i<nTables; ++i) {
226 if(!theT[i]) { continue; }
227 //G4cout << "## PreparePhysTable " << i << "." << G4endl;
228 G4PhysicsTable* table = theHandler->MakeTable(i);
229 //G4cout << " make table " << table << G4endl;
230 for(std::size_t j=0; j<numOfCouples; ++j) {
231 vec = (*table)[j];
232 if (bld->GetFlag(j) && nullptr == vec) {
233 //G4cout <<"InitialiseProcess iTable="<<i<<" jCouple="<< 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) {
261 thePhotoElectric->SetEmMasterProcess(theHandler->GetMasterProcess(0));
262 baseMat = theHandler->GetMasterProcess(0)->UseBaseMaterial();
263 }
264 thePhotoElectric->BuildPhysicsTable(part);
265
266 if(!isTheMaster) {
267 theCompton->SetEmMasterProcess(theHandler->GetMasterProcess(1));
268 }
269 theCompton->BuildPhysicsTable(part);
270
271 if(!isTheMaster) {
272 theConversionEE->SetEmMasterProcess(theHandler->GetMasterProcess(2));
273 }
274 theConversionEE->BuildPhysicsTable(part);
275
276 if(theRayleigh != nullptr) {
277 if(!isTheMaster) {
278 theRayleigh->SetEmMasterProcess(theHandler->GetMasterProcess(3));
279 }
280 theRayleigh->BuildPhysicsTable(part);
281 }
282 if(theGammaNuclear != nullptr) { theGammaNuclear->BuildPhysicsTable(part); }
283 if(theConversionMM != nullptr) { theConversionMM->BuildPhysicsTable(part); }
284
285 if(isTheMaster) {
286 const G4ProductionCutsTable* theCoupleTable=
288 G4int numOfCouples = (G4int)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(G4int i=0; i<numOfCouples; ++i) {
302
303 if (bld->GetFlag(i)) {
304 G4int idx = (!baseMat) ? i : DensityIndex(i);
305 const G4MaterialCutsCouple* couple =
306 theCoupleTable->GetMaterialCutsCouple(i);
307 const G4Material* material = couple->GetMaterial();
308
309 // energy interval 0
310 std::size_t nn = (*(tables[0]))[idx]->GetVectorLength();
311 if(1 < verboseLevel) {
312 G4cout << "======= Zone 0 ======= N= " << nn
313 << " for " << material->GetName() << G4endl;
314 }
315 for(std::size_t j=0; j<nn; ++j) {
316 G4double e = (*(tables[0]))[idx]->Energy(j);
317 G4double loge = G4Log(e);
318 sigComp = theCompton->GetLambda(e, couple, loge);
319 sigR = (nullptr != theRayleigh) ?
320 theRayleigh->GetLambda(e, couple, loge) : 0.0;
321 G4double sum = sigComp + sigR;
322 if(1 < verboseLevel) {
323 G4cout << j << ". E= " << e << " xs= " << sum
324 << " compt= " << sigComp << " Rayl= " << sigR << G4endl;
325 }
326 (*(tables[0]))[idx]->PutValue(j, sum);
327 if(theT[1]) {
328 val = sigR/sum;
329 (*(tables[1]))[idx]->PutValue(j, val);
330 }
331 }
332
333 // energy interval 1
334 nn = (*(tables[2]))[idx]->GetVectorLength();
335 if(1 < verboseLevel) {
336 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
337 }
338 for(std::size_t j=0; j<nn; ++j) {
339 G4double e = (*(tables[2]))[idx]->Energy(j);
340 G4double loge = G4Log(e);
341 sigComp = theCompton->GetLambda(e, couple, loge);
342 sigR = (nullptr != theRayleigh) ?
343 theRayleigh->GetLambda(e, couple, loge) : 0.0;
344 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
345 G4double sum = sigComp + sigR + sigPE;
346 if(1 < verboseLevel) {
347 G4cout << j << ". E= " << e << " xs= " << sum
348 << " compt= " << sigComp << " conv= " << sigConv
349 << " PE= " << sigPE << " Rayl= " << sigR
350 << " GN= " << sigN << G4endl;
351 }
352 (*(tables[2]))[idx]->PutValue(j, sum);
353
354 val = sigPE/sum;
355 (*(tables[3]))[idx]->PutValue(j, val);
356
357 val = (sigR > 0.0) ? (sigComp + sigPE)/sum : 1.0;
358 (*(tables[4]))[idx]->PutValue(j, val);
359 }
360
361 // energy interval 2
362 nn = (*(tables[6]))[idx]->GetVectorLength();
363 if(1 < verboseLevel) {
364 G4cout << "======= Zone 2 ======= N= " << nn << G4endl;
365 }
366 for(std::size_t j=0; j<nn; ++j) {
367 G4double e = (*(tables[6]))[idx]->Energy(j);
368 G4double loge = G4Log(e);
369 sigComp = theCompton->GetLambda(e, couple, loge);
370 sigConv = theConversionEE->GetLambda(e, couple, loge);
371 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
372 sigN = 0.0;
373 if(nullptr != gn) {
374 dynParticle->SetKineticEnergy(e);
375 sigN = gn->ComputeCrossSection(dynParticle, material);
376 }
377 G4double sum = sigComp + sigConv + sigPE + sigN;
378 if(1 < verboseLevel) {
379 G4cout << j << ". E= " << e << " xs= " << sum
380 << " compt= " << sigComp << " conv= " << sigConv
381 << " PE= " << sigPE
382 << " GN= " << sigN << G4endl;
383 }
384 (*(tables[6]))[idx]->PutValue(j, sum);
385
386 val = sigConv/sum;
387 (*(tables[7]))[idx]->PutValue(j, val);
388
389 val = (sigConv + sigComp)/sum;
390 (*(tables[8]))[idx]->PutValue(j, val);
391
392 val = (sigN > 0.0) ? (sigConv + sigComp + sigPE)/sum : 1.0;
393 (*(tables[9]))[idx]->PutValue(j, val);
394 }
395
396 // energy interval 3
397 nn = (*(tables[10]))[idx]->GetVectorLength();
398 if(1 < verboseLevel) {
399 G4cout << "======= Zone 3 ======= N= " << nn
400 << " for " << material->GetName() << G4endl;
401 }
402 for(std::size_t j=0; j<nn; ++j) {
403 G4double e = (*(tables[10]))[idx]->Energy(j);
404 G4double loge = G4Log(e);
405 sigComp = theCompton->GetLambda(e, couple, loge);
406 sigConv = theConversionEE->GetLambda(e, couple, loge);
407 sigPE = thePhotoElectric->GetLambda(e, couple, loge);
408 sigN = 0.0;
409 if(nullptr != gn) {
410 dynParticle->SetKineticEnergy(e);
411 sigN = gn->ComputeCrossSection(dynParticle, material);
412 }
413 sigM = 0.0;
414 if(nullptr != theConversionMM) {
415 val = theConversionMM->ComputeMeanFreePath(e, material);
416 sigM = (val < DBL_MAX) ? 1./val : 0.0;
417 }
418 G4double sum = sigComp + sigConv + sigPE + sigN + sigM;
419 if(1 < verboseLevel) {
420 G4cout << j << ". E= " << e << " xs= " << sum
421 << " compt= " << sigComp << " conv= " << sigConv
422 << " PE= " << sigPE
423 << " GN= " << sigN << G4endl;
424 }
425 (*(tables[10]))[idx]->PutValue(j, sum);
426
427 val = (sigComp + sigPE + sigN + sigM)/sum;
428 (*(tables[11]))[idx]->PutValue(j, val);
429
430 val = (sigPE + sigN + sigM)/sum;
431 (*(tables[12]))[idx]->PutValue(j, val);
432
433 val = (sigN + sigM)/sum;
434 (*(tables[13]))[idx]->PutValue(j, val);
435
436 val = sigN/sum;
437 (*(tables[14]))[idx]->PutValue(j, val);
438 }
439 for(std::size_t k=0; k<nTables; ++k) {
440 if(theT[k] && (k <= 1 || k >= 10)) {
441 //G4cout <<"BuildPhysTable spline iTable="<<k<<" jCouple="<< idx << G4endl;
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
475 G4double energy = track.GetKineticEnergy();
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) {
492 preStepKinEnergy = energy;
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 */
546 if(preStepKinEnergy < minPEEnergy) {
547 cross = ComputeGeneralLambda(0, 0);
548 //G4cout << "XS1: " << cross << G4endl;
549 peLambda = thePhotoElectric->GetLambda(preStepKinEnergy, currentCouple, preStepLogE);
550 cross += peLambda;
551 //G4cout << "XS2: " << peLambda << 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= " << preStepKinEnergy << 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) {
590 SelectEmProcess(step, thePhotoElectric);
591 } else {
592 if(theT[1] && preStepLambda*q < (preStepLambda - peLambda)*GetProbability(1) + peLambda) {
593 SelectEmProcess(step, theRayleigh);
594 } else {
595 SelectEmProcess(step, theCompton);
596 }
597 }
598 break;
599
600 case 1:
601 if(q <= GetProbability(3)) {
602 SelectEmProcess(step, thePhotoElectric);
603 } else if(q <= GetProbability(4)) {
604 SelectEmProcess(step, theCompton);
605 } else if(theRayleigh) {
606 SelectEmProcess(step, theRayleigh);
607 } else {
608 SelectEmProcess(step, thePhotoElectric);
609 }
610 break;
611
612 case 2:
613 if(q <= GetProbability(7)) {
614 SelectEmProcess(step, theConversionEE);
615 } else if(q <= GetProbability(8)) {
616 SelectEmProcess(step, theCompton);
617 } else if(q <= GetProbability(9)) {
618 SelectEmProcess(step, thePhotoElectric);
619 } else if(theGammaNuclear) {
620 SelectHadProcess(track, step, theGammaNuclear);
621 } else {
622 SelectEmProcess(step, theConversionEE);
623 }
624 break;
625
626 case 3:
627 if(q + GetProbability(11) <= 1.0) {
628 SelectEmProcess(step, theConversionEE);
629 } else if(q + GetProbability(12) <= 1.0) {
630 SelectEmProcess(step, theCompton);
631 } else if(q + GetProbability(13) <= 1.0) {
632 SelectEmProcess(step, thePhotoElectric);
633 } else if(theGammaNuclear && q + GetProbability(14) <= 1.0) {
634 SelectHadProcess(track, step, theGammaNuclear);
635 } else if(theConversionMM) {
636 SelectedProcess(step, theConversionMM);
637 } else {
638 SelectEmProcess(step, theConversionEE);
639 }
640 break;
641 }
642 // sample secondaries
643 if(selectedProc != nullptr) {
644 return selectedProc->PostStepDoIt(track, step);
645 }
646 // no interaction - exception case
648 return &fParticleChange;
649}
650
651//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
652
654 const G4Step& step, G4HadronicProcess* proc)
655{
656 SelectedProcess(step, proc);
659}
660
661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662
664 const G4String& directory,
665 G4bool ascii)
666{
667 G4bool yes = true;
668 if(!isTheMaster) { return yes; }
669 if(!thePhotoElectric->StorePhysicsTable(part, directory, ascii))
670 { yes = false; }
671 if(!theCompton->StorePhysicsTable(part, directory, ascii))
672 { yes = false; }
673 if(!theConversionEE->StorePhysicsTable(part, directory, ascii))
674 { yes = false; }
675 if(theRayleigh != nullptr &&
676 !theRayleigh->StorePhysicsTable(part, directory, ascii))
677 { yes = false; }
678
679 for(std::size_t i=0; i<nTables; ++i) {
680 if(theT[i]) {
681 G4String nam = (0==i || 2==i || 6==i || 10==i)
682 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
683 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
684 if(!theHandler->StorePhysicsTable(i, part, fnam, ascii)) { yes = false; }
685 }
686 }
687 return yes;
688}
689
690//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
691
692G4bool
694 const G4String& directory,
695 G4bool ascii)
696{
697 if(1 < verboseLevel) {
698 G4cout << "G4GammaGeneralProcess::RetrievePhysicsTable() for "
699 << part->GetParticleName() << " and process "
700 << GetProcessName() << G4endl;
701 }
702 G4bool yes = true;
703 if(!thePhotoElectric->RetrievePhysicsTable(part, directory, ascii))
704 { yes = false; }
705 if(!theCompton->RetrievePhysicsTable(part, directory, ascii))
706 { yes = false; }
707 if(!theConversionEE->RetrievePhysicsTable(part, directory, ascii))
708 { yes = false; }
709 if(theRayleigh != nullptr &&
710 !theRayleigh->RetrievePhysicsTable(part, directory, ascii))
711 { yes = false; }
712
713 for(std::size_t i=0; i<nTables; ++i) {
714 if(theT[i]) {
715 G4String nam = (0==i || 2==i || 6==i || 10==i)
716 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
717 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
718 G4bool spline = (i <= 1 || i >= 10);
719 if(!theHandler->RetrievePhysicsTable(i, part, fnam, ascii, spline))
720 { yes = false; }
721 }
722 }
723 if(yes) {
724 }
725 return yes;
726}
727
728//....Ooooo0ooooo ........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
729
731 G4double,
733{
735 return MeanFreePath(track);
736}
737
738//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
739
740void G4GammaGeneralProcess::ProcessDescription(std::ostream& out) const
741{
742 thePhotoElectric->ProcessDescription(out);
743 theCompton->ProcessDescription(out);
744 theConversionEE->ProcessDescription(out);
745 if(theRayleigh) { theRayleigh->ProcessDescription(out); }
747 if(theConversionMM) { theConversionMM->ProcessDescription(out); }
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
751
753{
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
761{
764}
765
766//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
767
769{
770 G4VEmProcess* proc = nullptr;
771 if(name == thePhotoElectric->GetProcessName()) {
772 proc = thePhotoElectric;
773 } else if(name == theCompton->GetProcessName()) {
774 proc = theCompton;
775 } else if(name == theConversionEE->GetProcessName()) {
776 proc = theConversionEE;
777 } else if(theRayleigh != nullptr && name == theRayleigh->GetProcessName()) {
778 proc = theRayleigh;
779 }
780 return proc;
781}
782
783//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
784
786{
787 return selectedProc;
788}
789
790//....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:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4double G4Log(G4double x)
Definition: G4Log.hh:227
@ 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: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 *)
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
const G4VProcess * GetCreatorProcess() const override
void InitialiseProcess(const G4ParticleDefinition *) override
void SelectEmProcess(const G4Step &, G4VEmProcess *)
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: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
const G4String & GetName() const
Definition: G4Material.hh:172
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 GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
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
G4ParticleChangeForGamma fParticleChange
void SetParticle(const G4ParticleDefinition *p)
void PreparePhysicsTable(const G4ParticleDefinition &) override
const G4Material * currentMaterial
G4double currentInteractionLength
Definition: G4VProcess.hh:339
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:416
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:181
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4int verboseLevel
Definition: G4VProcess.hh:360
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:188
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
Definition: DoubConv.h:17
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62