Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmModelManager Class Reference

#include <G4EmModelManager.hh>

Public Member Functions

 G4EmModelManager ()
 
 ~G4EmModelManager ()
 
void Clear ()
 
const G4DataVectorInitialise (const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4int)
 
void FillDEDXVector (G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
 
void FillLambdaVector (G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
 
G4VEmModelGetModel (G4int, G4bool ver=false)
 
void AddEmModel (G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
void DumpModelList (G4int verb)
 
G4VEmModelSelectModel (G4double &energy, size_t &index)
 
const G4DataVectorCuts () const
 
const G4DataVectorSubCutoff () const
 
void SetFluoFlag (G4bool val)
 
G4int NumberOfModels () const
 

Detailed Description

Definition at line 142 of file G4EmModelManager.hh.

Constructor & Destructor Documentation

◆ G4EmModelManager()

G4EmModelManager::G4EmModelManager ( )

Definition at line 124 of file G4EmModelManager.cc.

124 :
125 nEmModels(0),
126 nRegions(0),
127 particle(0),
128 verboseLevel(0)
129{
130 maxSubCutInRange = 0.7*mm;
131 models.reserve(4);
132 flucModels.reserve(4);
133 regions.reserve(4);
134 orderOfModels.reserve(4);
135 isUsed.reserve(4);
136 severalModels = true;
137 fluoFlag = false;
138 currRegionModel = 0;
139 currModel = 0;
140 theCuts = 0;
141 theSubCuts = 0;
142}

◆ ~G4EmModelManager()

G4EmModelManager::~G4EmModelManager ( )

Definition at line 146 of file G4EmModelManager.cc.

147{
148 verboseLevel = 0; // no verbosity at destruction
149 Clear();
150 delete theSubCuts;
151}

Member Function Documentation

◆ AddEmModel()

void G4EmModelManager::AddEmModel ( G4int  num,
G4VEmModel p,
G4VEmFluctuationModel fm,
const G4Region r 
)

Definition at line 171 of file G4EmModelManager.cc.

173{
174 if(!p) {
175 G4cout << "G4EmModelManager::AddEmModel WARNING: no model defined."
176 << G4endl;
177 return;
178 }
179 models.push_back(p);
180 flucModels.push_back(fm);
181 regions.push_back(r);
182 orderOfModels.push_back(num);
183 isUsed.push_back(0);
184 p->DefineForRegion(r);
185 ++nEmModels;
186}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:249

Referenced by G4VMultipleScattering::AddEmModel(), G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), and G4AdjointBremsstrahlungModel::G4AdjointBremsstrahlungModel().

◆ Clear()

void G4EmModelManager::Clear ( )

Definition at line 155 of file G4EmModelManager.cc.

156{
157 if(1 < verboseLevel) {
158 G4cout << "G4EmModelManager::Clear()" << G4endl;
159 }
160 size_t n = setOfRegionModels.size();
161 if(n > 0) {
162 for(size_t i=0; i<n; ++i) {
163 delete setOfRegionModels[i];
164 setOfRegionModels[i] = 0;
165 }
166 }
167}

Referenced by Initialise(), and ~G4EmModelManager().

◆ Cuts()

const G4DataVector * G4EmModelManager::Cuts ( ) const
inline

Definition at line 243 of file G4EmModelManager.hh.

244{
245 return theCuts;
246}

◆ DumpModelList()

void G4EmModelManager::DumpModelList ( G4int  verb)

Definition at line 663 of file G4EmModelManager.cc.

664{
665 if(verb == 0) { return; }
666 for(G4int i=0; i<nRegions; ++i) {
667 G4RegionModels* r = setOfRegionModels[i];
668 const G4Region* reg = r->Region();
669 G4int n = r->NumberOfModels();
670 if(n > 0) {
671 G4cout << " ===== EM models for the G4Region " << reg->GetName()
672 << " ======" << G4endl;;
673 for(G4int j=0; j<n; ++j) {
674 G4VEmModel* model = models[r->ModelIndex(j)];
675 G4double emin =
676 std::max(r->LowEdgeEnergy(j),model->LowEnergyActivationLimit());
677 G4double emax =
678 std::min(r->LowEdgeEnergy(j+1),model->HighEnergyActivationLimit());
679 G4cout << std::setw(20);
680 G4cout << model->GetName() << " : Emin= "
681 << std::setw(8) << G4BestUnit(emin,"Energy")
682 << " Emax= "
683 << std::setw(8) << G4BestUnit(emax,"Energy");
684 G4PhysicsTable* table = model->GetCrossSectionTable();
685 if(table) {
686 size_t kk = table->size();
687 for(size_t k=0; k<kk; ++k) {
688 G4PhysicsVector* v = (*table)[k];
689 if(v) {
690 G4int nn = v->GetVectorLength() - 1;
691 G4cout << " Table with " << nn << " bins Emin= "
692 << std::setw(6) << G4BestUnit(v->Energy(0),"Energy")
693 << " Emax= "
694 << std::setw(6) << G4BestUnit(v->Energy(nn),"Energy");
695 break;
696 }
697 }
698 }
700 if(an) { G4cout << " " << an->GetName(); }
701 if(fluoFlag && model->DeexcitationFlag()) { G4cout << " FluoActive"; }
702 G4cout << G4endl;
703 }
704 }
705 if(1 == nEmModels) { break; }
706 }
707}
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
size_t GetVectorLength() const
G4double Energy(size_t index) const
const G4String & GetName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:508
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:536
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:571
const G4String & GetName() const
Definition: G4VEmModel.hh:655
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:660
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:543

Referenced by G4VMultipleScattering::BuildPhysicsTable(), G4VEmProcess::PrintInfoDefinition(), G4VEnergyLossProcess::PrintInfoDefinition(), and G4VMultipleScattering::PrintInfoDefinition().

◆ FillDEDXVector()

void G4EmModelManager::FillDEDXVector ( G4PhysicsVector aVector,
const G4MaterialCutsCouple couple,
G4EmTableType  t = fRestricted 
)

Definition at line 510 of file G4EmModelManager.cc.

513{
514 size_t i = couple->GetIndex();
515 G4double cut = (*theCuts)[i];
516 G4double emin = 0.0;
517
518 if(fTotal == tType) { cut = DBL_MAX; }
519 else if(fSubRestricted == tType) {
520 emin = cut;
521 if(theSubCuts) { emin = (*theSubCuts)[i]; }
522 }
523
524 if(1 < verboseLevel) {
525 G4cout << "G4EmModelManager::FillDEDXVector() for "
526 << couple->GetMaterial()->GetName()
527 << " cut(MeV)= " << cut
528 << " emin(MeV)= " << emin
529 << " Type " << tType
530 << " for " << particle->GetParticleName()
531 << G4endl;
532 }
533
534 G4int reg = 0;
535 if(nRegions > 1 && nEmModels > 1) { reg = idxOfRegionModels[i]; }
536 const G4RegionModels* regModels = setOfRegionModels[reg];
537 G4int nmod = regModels->NumberOfModels();
538
539 // Calculate energy losses vector
540
541 //G4cout << "nmod= " << nmod << G4endl;
542 size_t totBinsLoss = aVector->GetVectorLength();
543 G4double del = 0.0;
544 G4int k0 = 0;
545
546 for(size_t j=0; j<totBinsLoss; ++j) {
547
548 G4double e = aVector->Energy(j);
549
550 // Choose a model of energy losses
551 G4int k = 0;
552 if (nmod > 1) {
553 k = nmod;
554 do {--k;} while (k>0 && e <= regModels->LowEdgeEnergy(k));
555 //G4cout << "k= " << k << G4endl;
556 if(k > 0 && k != k0) {
557 k0 = k;
558 G4double elow = regModels->LowEdgeEnergy(k);
559 G4double dedx1 = ComputeDEDX(models[regModels->ModelIndex(k-1)],
560 couple,elow,cut,emin);
561 G4double dedx2 = ComputeDEDX(models[regModels->ModelIndex(k)],
562 couple,elow,cut,emin);
563 del = 0.0;
564 if(dedx2 > 0.0) { del = (dedx1/dedx2 - 1.0)*elow; }
565 //G4cout << "elow= " << elow
566 // << " dedx1= " << dedx1 << " dedx2= " << dedx2 << G4endl;
567 }
568 }
569 G4double dedx =
570 ComputeDEDX(models[regModels->ModelIndex(k)],couple,e,cut,emin);
571 dedx *= (1.0 + del/e);
572
573 if(2 < verboseLevel) {
574 G4cout << "Material= " << couple->GetMaterial()->GetName()
575 << " E(MeV)= " << e/MeV
576 << " dEdx(MeV/mm)= " << dedx*mm/MeV
577 << " del= " << del*mm/MeV<< " k= " << k
578 << " modelIdx= " << regModels->ModelIndex(k)
579 << G4endl;
580 }
581 if(dedx < 0.0) { dedx = 0.0; }
582 aVector->PutValue(j, dedx);
583 }
584}
@ fSubRestricted
@ fTotal
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
const G4String & GetParticleName() const
void PutValue(size_t index, G4double theValue)
#define DBL_MAX
Definition: templates.hh:83

Referenced by G4VEnergyLossProcess::BuildDEDXTable().

◆ FillLambdaVector()

void G4EmModelManager::FillLambdaVector ( G4PhysicsVector aVector,
const G4MaterialCutsCouple couple,
G4bool  startFromNull = true,
G4EmTableType  t = fRestricted 
)

Definition at line 588 of file G4EmModelManager.cc.

592{
593 size_t i = couple->GetIndex();
594 G4double cut = (*theCuts)[i];
595 G4double tmax = DBL_MAX;
596 if (fSubRestricted == tType) {
597 tmax = cut;
598 if(theSubCuts) { cut = (*theSubCuts)[i]; }
599 }
600
601 G4int reg = 0;
602 if(nRegions > 1 && nEmModels > 1) { reg = idxOfRegionModels[i]; }
603 const G4RegionModels* regModels = setOfRegionModels[reg];
604 G4int nmod = regModels->NumberOfModels();
605 if(1 < verboseLevel) {
606 G4cout << "G4EmModelManager::FillLambdaVector() for "
607 << particle->GetParticleName()
608 << " in " << couple->GetMaterial()->GetName()
609 << " Ecut(MeV)= " << cut
610 << " Emax(MeV)= " << tmax
611 << " Type " << tType
612 << " nmod= " << nmod
613 << G4endl;
614 }
615
616 // Calculate lambda vector
617 size_t totBinsLambda = aVector->GetVectorLength();
618 G4double del = 0.0;
619 G4int k0 = 0;
620 G4int k = 0;
621 G4VEmModel* mod = models[regModels->ModelIndex(0)];
622 for(size_t j=0; j<totBinsLambda; ++j) {
623
624 G4double e = aVector->Energy(j);
625
626 // Choose a model
627 if (nmod > 1) {
628 k = nmod;
629 do {--k;} while (k>0 && e <= regModels->LowEdgeEnergy(k));
630 if(k > 0 && k != k0) {
631 k0 = k;
632 G4double elow = regModels->LowEdgeEnergy(k);
633 G4VEmModel* mod1 = models[regModels->ModelIndex(k-1)];
634 G4double xs1 = mod1->CrossSection(couple,particle,elow,cut,tmax);
635 mod = models[regModels->ModelIndex(k)];
636 G4double xs2 = mod->CrossSection(couple,particle,elow,cut,tmax);
637 del = 0.0;
638 if(xs2 > 0.0) { del = (xs1/xs2 - 1.0)*elow; }
639 //G4cout << "New model k=" << k << " E(MeV)= " << e/MeV
640 // << " Elow(MeV)= " << elow/MeV << " del= " << del << G4endl;
641 }
642 }
643 G4double cross = mod->CrossSection(couple,particle,e,cut,tmax);
644 cross *= (1.0 + del/e);
645 if(fIsCrossSectionPrim == tType) { cross *= e; }
646
647 if(j==0 && startFromNull) { cross = 0.0; }
648
649 if(2 < verboseLevel) {
650 G4cout << "FillLambdaVector: " << j << ". e(MeV)= " << e/MeV
651 << " cross(1/mm)= " << cross*mm
652 << " del= " << del*mm << " k= " << k
653 << " modelIdx= " << regModels->ModelIndex(k)
654 << G4endl;
655 }
656 if(cross < 0.0) { cross = 0.0; }
657 aVector->PutValue(j, cross);
658 }
659}
@ fIsCrossSectionPrim
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:418

Referenced by G4VEnergyLossProcess::BuildLambdaTable().

◆ GetModel()

G4VEmModel * G4EmModelManager::GetModel ( G4int  i,
G4bool  ver = false 
)

Definition at line 209 of file G4EmModelManager.cc.

210{
211 G4VEmModel* model = 0;
212 if(i >= 0 && i < nEmModels) { model = models[i]; }
213 else if(verboseLevel > 0 && ver) {
214 G4cout << "G4EmModelManager::GetModel WARNING: "
215 << "index " << i << " is wrong Nmodels= "
216 << nEmModels;
217 if(particle) G4cout << " for " << particle->GetParticleName();
218 G4cout<< G4endl;
219 }
220 return model;
221}

Referenced by G4VEmProcess::GetModelByIndex(), G4VEnergyLossProcess::GetModelByIndex(), G4VMultipleScattering::GetModelByIndex(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4VMultipleScattering::SetIonisation(), G4VMultipleScattering::StartTracking(), and G4VMultipleScattering::StorePhysicsTable().

◆ Initialise()

const G4DataVector * G4EmModelManager::Initialise ( const G4ParticleDefinition p,
const G4ParticleDefinition secondaryParticle,
G4double  minSubRange,
G4int  val 
)

Definition at line 226 of file G4EmModelManager.cc.

230{
231 verboseLevel = val;
232 G4String partname = p->GetParticleName();
233 if(1 < verboseLevel) {
234 G4cout << "G4EmModelManager::Initialise() for "
235 << partname << G4endl;
236 }
237 // Are models defined?
238 if(nEmModels < 1) {
240 ed << "No models found out for " << p->GetParticleName()
241 << " !" << G4endl;
242 G4Exception("G4EmModelManager::Initialise","em0002",
243 FatalException, ed);
244 }
245
246 particle = p;
247 Clear(); // needed if run is not first
248
250 const G4Region* world =
251 regionStore->GetRegion("DefaultRegionForTheWorld", false);
252
253 // Identify the list of regions with different set of models
254 nRegions = 1;
255 std::vector<const G4Region*> setr;
256 setr.push_back(world);
257 G4bool isWorld = false;
258
259 for (G4int ii=0; ii<nEmModels; ++ii) {
260 const G4Region* r = regions[ii];
261 if ( r == 0 || r == world) {
262 isWorld = true;
263 regions[ii] = world;
264 } else {
265 G4bool newRegion = true;
266 if (nRegions>1) {
267 for (G4int j=1; j<nRegions; ++j) {
268 if ( r == setr[j] ) newRegion = false;
269 }
270 }
271 if (newRegion) {
272 setr.push_back(r);
273 nRegions++;
274 }
275 }
276 }
277 // Are models defined?
278 if(!isWorld) {
280 ed << "No models defined for the World volume for " << p->GetParticleName()
281 << " !" << G4endl;
282 G4Exception("G4EmModelManager::Initialise","em0002",
283 FatalException, ed);
284 }
285
286 G4ProductionCutsTable* theCoupleTable=
288 size_t numOfCouples = theCoupleTable->GetTableSize();
289
290 // prepare vectors, shortcut for the case of only 1 model
291 if(nRegions > 1 && nEmModels > 1) {
292 if(numOfCouples > idxOfRegionModels.size()) {
293 idxOfRegionModels.resize(numOfCouples);
294 }
295 }
296 size_t nr = 1;
297 if(nEmModels > 1) { nr = nRegions; }
298 if(nr > setOfRegionModels.size()) { setOfRegionModels.resize(nr); }
299
300 std::vector<G4int> modelAtRegion(nEmModels);
301 std::vector<G4int> modelOrd(nEmModels);
302 G4DataVector eLow(nEmModels+1);
303 G4DataVector eHigh(nEmModels);
304
305 // Order models for regions
306 for (G4int reg=0; reg<nRegions; ++reg) {
307 const G4Region* region = setr[reg];
308 G4int n = 0;
309
310 for (G4int ii=0; ii<nEmModels; ++ii) {
311
312 G4VEmModel* model = models[ii];
313 if ( region == regions[ii] ) {
314
315 G4double tmin = model->LowEnergyLimit();
316 G4double tmax = model->HighEnergyLimit();
317 G4int ord = orderOfModels[ii];
318 G4bool push = true;
319 G4bool insert = false;
320 G4int idx = n;
321
322 if(1 < verboseLevel) {
323 G4cout << "Model #" << ii
324 << " <" << model->GetName() << "> for region <";
325 if (region) G4cout << region->GetName();
326 G4cout << "> "
327 << " tmin(MeV)= " << tmin/MeV
328 << "; tmax(MeV)= " << tmax/MeV
329 << "; order= " << ord
330 << G4endl;
331 }
332
333 if(n > 0) {
334
335 // extend energy range to previous models
336 tmin = std::min(tmin, eHigh[n-1]);
337 tmax = std::max(tmax, eLow[0]);
338 //G4cout << "tmin= " << tmin << " tmax= "
339 // << tmax << " ord= " << ord <<G4endl;
340 // empty energy range
341 if( tmax - tmin <= eV) push = false;
342 // low-energy model
343 else if (tmax == eLow[0]) {
344 push = false;
345 insert = true;
346 idx = 0;
347 // resolve intersections
348 } else if(tmin < eHigh[n-1]) {
349 // compare order
350 for(G4int k=0; k<n; ++k) {
351 // new model has lower application
352 if(ord >= modelOrd[k]) {
353 if(tmin < eHigh[k] && tmin >= eLow[k]) tmin = eHigh[k];
354 if(tmax <= eHigh[k] && tmax > eLow[k]) tmax = eLow[k];
355 if(tmax > eHigh[k] && tmin < eLow[k]) {
356 if(tmax - eHigh[k] > eLow[k] - tmin) tmin = eHigh[k];
357 else tmax = eLow[k];
358 }
359 if( tmax - tmin <= eV) {
360 push = false;
361 break;
362 }
363 }
364 }
365 //G4cout << "tmin= " << tmin << " tmax= "
366 // << tmax << " push= " << push << " idx= " << idx <<G4endl;
367 if(push) {
368 if (tmax == eLow[0]) {
369 push = false;
370 insert = true;
371 idx = 0;
372 // continue resolve intersections
373 } else if(tmin < eHigh[n-1]) {
374 // last energy interval
375 if(tmin > eLow[n-1] && tmax >= eHigh[n-1]) {
376 eHigh[n-1] = tmin;
377 // first energy interval
378 } else if(tmin <= eLow[0] && tmax < eHigh[0]) {
379 eLow[0] = tmax;
380 push = false;
381 insert = true;
382 idx = 0;
383 } else {
384 // find energy interval to replace
385 for(G4int k=0; k<n; ++k) {
386 if(tmin <= eLow[k] && tmax >= eHigh[k]) {
387 push = false;
388 modelAtRegion[k] = ii;
389 modelOrd[k] = ord;
390 isUsed[ii] = 1;
391 }
392 }
393 }
394 }
395 }
396 }
397 }
398 if(insert) {
399 for(G4int k=n-1; k>=idx; --k) {
400 modelAtRegion[k+1] = modelAtRegion[k];
401 modelOrd[k+1] = modelOrd[k];
402 eLow[k+1] = eLow[k];
403 eHigh[k+1] = eHigh[k];
404 }
405 }
406 //G4cout << "push= " << push << " insert= " << insert
407 //<< " idx= " << idx <<G4endl;
408 if (push || insert) {
409 ++n;
410 modelAtRegion[idx] = ii;
411 modelOrd[idx] = ord;
412 eLow[idx] = tmin;
413 eHigh[idx] = tmax;
414 isUsed[ii] = 1;
415 }
416 }
417 }
418 eLow[0] = 0.0;
419 eLow[n] = eHigh[n-1];
420
421 if(1 < verboseLevel) {
422 G4cout << "New G4RegionModels set with " << n << " models for region <";
423 if (region) G4cout << region->GetName();
424 G4cout << "> Elow(MeV)= ";
425 for(G4int ii=0; ii<=n; ++ii) {G4cout << eLow[ii]/MeV << " ";}
426 G4cout << G4endl;
427 }
428 G4RegionModels* rm = new G4RegionModels(n, modelAtRegion, eLow, region);
429 setOfRegionModels[reg] = rm;
430 if(1 == nEmModels) { break; }
431 }
432
433 currRegionModel = setOfRegionModels[0];
434
435 // Access to materials and build cuts
436 size_t idx = 1;
437 if(secondaryParticle) {
438 if( secondaryParticle == G4Gamma::Gamma() ) { idx = 0; }
439 else if( secondaryParticle == G4Positron::Positron()) { idx = 2; }
440 }
441
442 //theCuts = theCoupleTable->GetEnergyCutsVector(idx);
443 theCuts = static_cast<const G4DataVector*>(theCoupleTable->GetEnergyCutsVector(idx));
444
445 if(minSubRange < 1.0) {
446 if( !theSubCuts ) { theSubCuts = new G4DataVector(); }
447 theSubCuts->resize(numOfCouples,DBL_MAX);
448 }
449 for(size_t i=0; i<numOfCouples; ++i) {
450
451 const G4MaterialCutsCouple* couple =
452 theCoupleTable->GetMaterialCutsCouple(i);
453 const G4Material* material = couple->GetMaterial();
454 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
455
456 G4int reg = 0;
457 if(nRegions > 1 && nEmModels > 1) {
458 reg = nRegions;
459 do {--reg;} while (reg>0 && pcuts != (setr[reg]->GetProductionCuts()));
460 idxOfRegionModels[i] = reg;
461 }
462 if(1 < verboseLevel) {
463 G4cout << "G4EmModelManager::Initialise() for "
464 << material->GetName()
465 << " indexOfCouple= " << i
466 << " indexOfRegion= " << reg
467 << G4endl;
468 }
469
470 G4double cut = (*theCuts)[i];
471 if(secondaryParticle) {
472
473 // compute subcut
474 if( cut < DBL_MAX && minSubRange < 1.0) {
475 G4double subcut = minSubRange*cut;
476 G4double rcut = std::min(minSubRange*pcuts->GetProductionCut(idx),
477 maxSubCutInRange);
478 G4double tcutmax =
479 theCoupleTable->ConvertRangeToEnergy(secondaryParticle,material,rcut);
480 if(tcutmax < subcut) { subcut = tcutmax; }
481 (*theSubCuts)[i] = subcut;
482 }
483 }
484 }
485
486 // initialize models
487 G4int nn = 0;
488 severalModels = true;
489 for(G4int jj=0; jj<nEmModels; ++jj) {
490 if(1 == isUsed[jj]) {
491 ++nn;
492 currModel = models[jj];
493 currModel->Initialise(particle, *theCuts);
494 if(flucModels[jj]) { flucModels[jj]->InitialiseMe(particle); }
495 }
496 }
497 if(1 == nn) { severalModels = false; }
498
499 if(1 < verboseLevel) {
500 G4cout << "G4EmModelManager for " << particle->GetParticleName()
501 << " is initialised; nRegions= " << nRegions
502 << G4endl;
503 }
504
505 return theCuts;
506}
@ FatalException
bool G4bool
Definition: G4Types.hh:67
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ProductionCuts * GetProductionCuts() const
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double ConvertRangeToEnergy(const G4ParticleDefinition *particle, const G4Material *material, G4double range)
G4double GetProductionCut(G4int index) const
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Referenced by G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecond(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), and G4VMultipleScattering::PreparePhysicsTable().

◆ NumberOfModels()

◆ SelectModel()

G4VEmModel * G4EmModelManager::SelectModel ( G4double energy,
size_t &  index 
)
inline

Definition at line 229 of file G4EmModelManager.hh.

231{
232 if(severalModels) {
233 if(nRegions > 1) {
234 currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
235 }
236 currModel = models[currRegionModel->SelectIndex(kinEnergy)];
237 }
238 return currModel;
239}

Referenced by G4VEmProcess::SelectModel(), G4VEnergyLossProcess::SelectModel(), G4VMultipleScattering::SelectModel(), G4VEnergyLossProcess::SelectModelForMaterial(), and G4VEmProcess::SelectModelForMaterial().

◆ SetFluoFlag()

void G4EmModelManager::SetFluoFlag ( G4bool  val)
inline

Definition at line 257 of file G4EmModelManager.hh.

258{
259 fluoFlag = val;
260}

Referenced by G4VEmProcess::PreparePhysicsTable().

◆ SubCutoff()

const G4DataVector * G4EmModelManager::SubCutoff ( ) const
inline

Definition at line 250 of file G4EmModelManager.hh.

251{
252 return theSubCuts;
253}

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ UpdateEmModel()

void G4EmModelManager::UpdateEmModel ( const G4String nam,
G4double  emin,
G4double  emax 
)

Definition at line 190 of file G4EmModelManager.cc.

192{
193 if (nEmModels > 0) {
194 for(G4int i=0; i<nEmModels; ++i) {
195 if(nam == models[i]->GetName()) {
196 models[i]->SetLowEnergyLimit(emin);
197 models[i]->SetHighEnergyLimit(emax);
198 break;
199 }
200 }
201 }
202 G4cout << "G4EmModelManager::UpdateEmModel WARNING: no model <"
203 << nam << "> is found out"
204 << G4endl;
205}

Referenced by G4VEmProcess::UpdateEmModel(), and G4VEnergyLossProcess::UpdateEmModel().


The documentation for this class was generated from the following files: