Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmModelManager.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: G4EmModelManager
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 07.05.2002
37//
38// Modifications: V.Ivanchenko
39//
40// Class Description:
41//
42// It is the unified energy loss process it calculates the continuous
43// energy loss for charged particles using a set of Energy Loss
44// models valid for different energy regions. There are a possibility
45// to create and access to dE/dx and range tables, or to calculate
46// that information on fly.
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
52#include "G4EmModelManager.hh"
53#include "G4SystemOfUnits.hh"
54#include "G4PhysicsTable.hh"
55#include "G4PhysicsVector.hh"
56#include "G4VMscModel.hh"
57
58#include "G4Step.hh"
60#include "G4PhysicsVector.hh"
63#include "G4RegionStore.hh"
64#include "G4Gamma.hh"
65#include "G4Electron.hh"
66#include "G4Positron.hh"
67#include "G4UnitsTable.hh"
68#include "G4DataVector.hh"
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
72G4RegionModels::G4RegionModels(G4int nMod, std::vector<G4int>& indx,
73 G4DataVector& lowE, const G4Region* reg)
74{
75 nModelsForRegion = nMod;
76 theListOfModelIndexes = new G4int [nModelsForRegion];
77 lowKineticEnergy = new G4double [nModelsForRegion+1];
78 for (G4int i=0; i<nModelsForRegion; ++i) {
79 theListOfModelIndexes[i] = indx[i];
80 lowKineticEnergy[i] = lowE[i];
81 }
82 lowKineticEnergy[nModelsForRegion] = lowE[nModelsForRegion];
83 theRegion = reg;
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
88G4RegionModels::~G4RegionModels()
89{
90 delete [] theListOfModelIndexes;
91 delete [] lowKineticEnergy;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
98{
99 models.reserve(4);
100 flucModels.reserve(4);
101 regions.reserve(4);
102 orderOfModels.reserve(4);
103 isUsed.reserve(4);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109{
110 verboseLevel = 0; // no verbosity at destruction
111 Clear();
112 delete theCutsNew;
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118{
119 if(1 < verboseLevel) {
120 G4cout << "G4EmModelManager::Clear()" << G4endl;
121 }
122 std::size_t n = setOfRegionModels.size();
123 for(std::size_t i=0; i<n; ++i) {
124 delete setOfRegionModels[i];
125 setOfRegionModels[i] = nullptr;
126 }
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130
132 G4VEmFluctuationModel* fm, const G4Region* r)
133{
134 if(nullptr == p) {
135 G4cout << "G4EmModelManager::AddEmModel WARNING: no model defined."
136 << G4endl;
137 return;
138 }
139 models.push_back(p);
140 flucModels.push_back(fm);
141 regions.push_back(r);
142 orderOfModels.push_back(num);
143 isUsed.push_back(0);
144 p->DefineForRegion(r);
145 ++nEmModels;
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
151{
152 G4VEmModel* model = nullptr;
153 if(idx >= 0 && idx < nEmModels) { model = models[idx]; }
154 else if(verboseLevel > 0 && ver) {
155 G4cout << "G4EmModelManager::GetModel WARNING: "
156 << "index " << idx << " is wrong Nmodels= "
157 << nEmModels;
158 if(nullptr != particle) {
159 G4cout << " for " << particle->GetParticleName();
160 }
161 G4cout<< G4endl;
162 }
163 return model;
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167
169{
170 G4RegionModels* rm = setOfRegionModels[idxOfRegionModels[idx]];
171 return (k < rm->NumberOfModels()) ? models[rm->ModelIndex(k)] : nullptr;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175
177{
178 G4RegionModels* rm = setOfRegionModels[idxOfRegionModels[idx]];
179 return rm->NumberOfModels();
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
184const G4DataVector*
186 const G4ParticleDefinition* secondaryParticle,
187 G4int verb)
188{
189 verboseLevel = verb;
190 if(1 < verboseLevel) {
191 G4cout << "G4EmModelManager::Initialise() for "
192 << p->GetParticleName() << " Nmodels= " << nEmModels << G4endl;
193 }
194 // Are models defined?
195 if(nEmModels < 1) {
197 ed << "No models found out for " << p->GetParticleName()
198 << " !";
199 G4Exception("G4EmModelManager::Initialise","em0002",
200 FatalException, ed);
201 }
202
203 particle = p;
204 Clear(); // needed if run is not first
205
207 const G4Region* world =
208 regionStore->GetRegion("DefaultRegionForTheWorld", false);
209
210 // Identify the list of regions with different set of models
211 nRegions = 1;
212 std::vector<const G4Region*> setr;
213 setr.push_back(world);
214 G4bool isWorld = false;
215
216 for (G4int ii=0; ii<nEmModels; ++ii) {
217 const G4Region* r = regions[ii];
218 if ( r == nullptr || r == world) {
219 isWorld = true;
220 regions[ii] = world;
221 } else {
222 G4bool newRegion = true;
223 if (nRegions>1) {
224 for (G4int j=1; j<nRegions; ++j) {
225 if ( r == setr[j] ) { newRegion = false; }
226 }
227 }
228 if (newRegion) {
229 setr.push_back(r);
230 ++nRegions;
231 }
232 }
233 }
234 // Are models defined?
235 if(!isWorld) {
237 ed << "No models defined for the World volume for "
238 << p->GetParticleName() << " !";
239 G4Exception("G4EmModelManager::Initialise","em0002",
240 FatalException, ed);
241 }
242
243 G4ProductionCutsTable* theCoupleTable=
245 std::size_t numOfCouples = theCoupleTable->GetTableSize();
246
247 // prepare vectors, shortcut for the case of only 1 model
248 // or only one region
249 if(nRegions > 1 && nEmModels > 1) {
250 idxOfRegionModels.resize(numOfCouples,0);
251 setOfRegionModels.resize((std::size_t)nRegions,nullptr);
252 } else {
253 idxOfRegionModels.resize(1,0);
254 setOfRegionModels.resize(1,nullptr);
255 }
256
257 std::vector<G4int> modelAtRegion(nEmModels);
258 std::vector<G4int> modelOrd(nEmModels);
259 G4DataVector eLow(nEmModels+1);
260 G4DataVector eHigh(nEmModels);
261
262 if(1 < verboseLevel) {
263 G4cout << " Nregions= " << nRegions
264 << " Nmodels= " << nEmModels << G4endl;
265 }
266
267 // Order models for regions
268 for (G4int reg=0; reg<nRegions; ++reg) {
269 const G4Region* region = setr[reg];
270 G4int n = 0;
271
272 for (G4int ii=0; ii<nEmModels; ++ii) {
273
274 G4VEmModel* model = models[ii];
275 if ( region == regions[ii] ) {
276
277 G4double tmin = model->LowEnergyLimit();
278 G4double tmax = model->HighEnergyLimit();
279 G4int ord = orderOfModels[ii];
280 G4bool push = true;
281 G4bool insert = false;
282 G4int idx = n;
283
284 if(1 < verboseLevel) {
285 G4cout << "Model #" << ii
286 << " <" << model->GetName() << "> for region <";
287 if (region) G4cout << region->GetName();
288 G4cout << "> "
289 << " tmin(MeV)= " << tmin/MeV
290 << "; tmax(MeV)= " << tmax/MeV
291 << "; order= " << ord
292 << "; tminAct= " << model->LowEnergyActivationLimit()/MeV
293 << "; tmaxAct= " << model->HighEnergyActivationLimit()/MeV
294 << G4endl;
295 }
296
297 static const G4double limitdelta = 0.01*eV;
298 if(n > 0) {
299
300 // extend energy range to previous models
301 tmin = std::min(tmin, eHigh[n-1]);
302 tmax = std::max(tmax, eLow[0]);
303 //G4cout << "tmin= " << tmin << " tmax= "
304 // << tmax << " ord= " << ord <<G4endl;
305 // empty energy range
306 if( tmax - tmin <= limitdelta) { push = false; }
307 // low-energy model
308 else if (tmax == eLow[0]) {
309 push = false;
310 insert = true;
311 idx = 0;
312 // resolve intersections
313 } else if(tmin < eHigh[n-1]) {
314 // compare order
315 for(G4int k=0; k<n; ++k) {
316 // new model has higher order parameter,
317 // so, its application area may be reduced
318 // to avoid intersections
319 if(ord >= modelOrd[k]) {
320 if(tmin < eHigh[k] && tmin >= eLow[k]) { tmin = eHigh[k]; }
321 if(tmax <= eHigh[k] && tmax > eLow[k]) { tmax = eLow[k]; }
322 if(tmax > eHigh[k] && tmin < eLow[k]) {
323 if(tmax - eHigh[k] > eLow[k] - tmin) { tmin = eHigh[k]; }
324 else { tmax = eLow[k]; }
325 }
326 if( tmax - tmin <= limitdelta) {
327 push = false;
328 break;
329 }
330 }
331 }
332 // this model has lower order parameter than possible
333 // other models, with which there may be intersections
334 // so, appliction area of such models may be reduced
335
336 // insert below the first model
337 if (tmax <= eLow[0]) {
338 push = false;
339 insert = true;
340 idx = 0;
341 // resolve intersections
342 } else if(tmin < eHigh[n-1]) {
343 // last energy interval
344 if(tmin > eLow[n-1] && tmax >= eHigh[n-1]) {
345 eHigh[n-1] = tmin;
346 // first energy interval
347 } else if(tmin <= eLow[0] && tmax < eHigh[0]) {
348 eLow[0] = tmax;
349 push = false;
350 insert = true;
351 idx = 0;
352 // loop over all models
353 } else {
354 for(G4int k=n-1; k>=0; --k) {
355 if(tmin <= eLow[k] && tmax >= eHigh[k]) {
356 // full overlap exclude previous model
357 isUsed[modelAtRegion[k]] = 0;
358 idx = k;
359 if(k < n-1) {
360 // shift upper models and change index
361 for(G4int kk=k; kk<n-1; ++kk) {
362 modelAtRegion[kk] = modelAtRegion[kk+1];
363 modelOrd[kk] = modelOrd[kk+1];
364 eLow[kk] = eLow[kk+1];
365 eHigh[kk] = eHigh[kk+1];
366 }
367 ++k;
368 }
369 --n;
370 } else {
371 // partially reduce previous model area
372 if(tmin <= eLow[k] && tmax > eLow[k]) {
373 eLow[k] = tmax;
374 idx = k;
375 insert = true;
376 push = false;
377 } else if(tmin < eHigh[k] && tmax >= eHigh[k]) {
378 eHigh[k] = tmin;
379 idx = k + 1;
380 if(idx < n) {
381 insert = true;
382 push = false;
383 }
384 } else if(tmin > eLow[k] && tmax < eHigh[k]) {
385 if(eHigh[k] - tmax > tmin - eLow[k]) {
386 eLow[k] = tmax;
387 idx = k;
388 insert = true;
389 push = false;
390 } else {
391 eHigh[k] = tmin;
392 idx = k + 1;
393 if(idx < n) {
394 insert = true;
395 push = false;
396 }
397 }
398 }
399 }
400 }
401 }
402 }
403 }
404 }
405 // provide space for the new model
406 if(insert) {
407 for(G4int k=n-1; k>=idx; --k) {
408 modelAtRegion[k+1] = modelAtRegion[k];
409 modelOrd[k+1] = modelOrd[k];
410 eLow[k+1] = eLow[k];
411 eHigh[k+1] = eHigh[k];
412 }
413 }
414 //G4cout << "push= " << push << " insert= " << insert
415 // << " idx= " << idx <<G4endl;
416 // the model is added
417 if (push || insert) {
418 ++n;
419 modelAtRegion[idx] = ii;
420 modelOrd[idx] = ord;
421 eLow[idx] = tmin;
422 eHigh[idx] = tmax;
423 isUsed[ii] = 1;
424 }
425 // exclude models with zero energy range
426 for(G4int k=n-1; k>=0; --k) {
427 if(eHigh[k] - eLow[k] <= limitdelta) {
428 isUsed[modelAtRegion[k]] = 0;
429 if(k < n-1) {
430 for(G4int kk=k; kk<n-1; ++kk) {
431 modelAtRegion[kk] = modelAtRegion[kk+1];
432 modelOrd[kk] = modelOrd[kk+1];
433 eLow[kk] = eLow[kk+1];
434 eHigh[kk] = eHigh[kk+1];
435 }
436 }
437 --n;
438 }
439 }
440 }
441 }
442 eLow[0] = 0.0;
443 eLow[n] = eHigh[n-1];
444
445 if(1 < verboseLevel) {
446 G4cout << "### New G4RegionModels set with " << n
447 << " models for region <";
448 if (region) { G4cout << region->GetName(); }
449 G4cout << "> Elow(MeV)= ";
450 for(G4int iii=0; iii<=n; ++iii) {G4cout << eLow[iii]/MeV << " ";}
451 G4cout << G4endl;
452 }
453 auto rm = new G4RegionModels(n, modelAtRegion, eLow, region);
454 setOfRegionModels[reg] = rm;
455 // shortcut
456 if(1 == nEmModels) { break; }
457 }
458
459 currRegionModel = setOfRegionModels[0];
460 currModel = models[0];
461
462 // Access to materials and build cuts
463 std::size_t idx = 1;
464 if(nullptr != secondaryParticle) {
465 if( secondaryParticle == G4Gamma::Gamma() ) { idx = 0; }
466 else if( secondaryParticle == G4Electron::Electron()) { idx = 1; }
467 else if( secondaryParticle == G4Positron::Positron()) { idx = 2; }
468 else { idx = 3; }
469 }
470
471 theCuts =
472 static_cast<const G4DataVector*>(theCoupleTable->GetEnergyCutsVector(idx));
473
474 // for the second run the check on cuts should be repeated
475 if(nullptr != theCutsNew) { *theCutsNew = *theCuts; }
476
477 // G4cout << "========Start define cuts" << G4endl;
478 // define cut values
479 for(std::size_t i=0; i<numOfCouples; ++i) {
480
481 const G4MaterialCutsCouple* couple =
482 theCoupleTable->GetMaterialCutsCouple((G4int)i);
483 const G4Material* material = couple->GetMaterial();
484 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
485
486 G4int reg = 0;
487 if(nRegions > 1 && nEmModels > 1) {
488 reg = nRegions;
489 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
490 do {--reg;} while (reg>0 && pcuts != (setr[reg]->GetProductionCuts()));
491 idxOfRegionModels[i] = reg;
492 }
493 if(1 < verboseLevel) {
494 G4cout << "G4EmModelManager::Initialise() for "
495 << material->GetName()
496 << " indexOfCouple= " << i
497 << " indexOfRegion= " << reg
498 << G4endl;
499 }
500
501 G4double cut = (*theCuts)[i];
502 if(nullptr != secondaryParticle) {
503
504 // note that idxOfRegionModels[] not always filled
505 G4int inn = 0;
506 G4int nnm = 1;
507 if(nRegions > 1 && nEmModels > 1) {
508 inn = idxOfRegionModels[i];
509 }
510 // check cuts and introduce upper limits
511 //G4cout << "idx= " << i << " cut(keV)= " << cut/keV << G4endl;
512 currRegionModel = setOfRegionModels[inn];
513 nnm = currRegionModel->NumberOfModels();
514
515 //G4cout << "idx= " << i << " Nmod= " << nnm << G4endl;
516
517 for(G4int jj=0; jj<nnm; ++jj) {
518 //G4cout << "jj= " << jj << " modidx= "
519 // << currRegionModel->ModelIndex(jj) << G4endl;
520 currModel = models[currRegionModel->ModelIndex(jj)];
521 G4double cutlim = currModel->MinEnergyCut(particle,couple);
522 if(cutlim > cut) {
523 if(nullptr == theCutsNew) { theCutsNew = new G4DataVector(*theCuts); }
524 (*theCutsNew)[i] = cutlim;
525 /*
526 G4cout << "### " << partname << " energy loss model in "
527 << material->GetName()
528 << " Cut was changed from " << cut/keV << " keV to "
529 << cutlim/keV << " keV " << " due to "
530 << currModel->GetName() << G4endl;
531 */
532 }
533 }
534 }
535 }
536 if(nullptr != theCutsNew) { theCuts = theCutsNew; }
537
538 // initialize models
539 G4int nn = 0;
540 severalModels = true;
541 for(G4int jj=0; jj<nEmModels; ++jj) {
542 if(1 == isUsed[jj]) {
543 ++nn;
544 currModel = models[jj];
545 currModel->Initialise(particle, *theCuts);
546 if(nullptr != flucModels[jj]) { flucModels[jj]->InitialiseMe(particle); }
547 }
548 }
549 if(1 == nn) { severalModels = false; }
550
551 if(1 < verboseLevel) {
552 G4cout << "G4EmModelManager for " << particle->GetParticleName()
553 << " is initialised; nRegions= " << nRegions
554 << " severalModels: " << severalModels
555 << G4endl;
556 }
557 return theCuts;
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
563 const G4MaterialCutsCouple* couple,
564 G4EmTableType tType)
565{
566 std::size_t i = couple->GetIndex();
567 G4double cut = (fTotal == tType) ? DBL_MAX : (*theCuts)[i];
568
569 if(1 < verboseLevel) {
570 G4cout << "G4EmModelManager::FillDEDXVector() for "
571 << couple->GetMaterial()->GetName()
572 << " cut(MeV)= " << cut
573 << " Type " << tType
574 << " for " << particle->GetParticleName()
575 << G4endl;
576 }
577
578 G4int reg = 0;
579 if(nRegions > 1 && nEmModels > 1) { reg = idxOfRegionModels[i]; }
580 const G4RegionModels* regModels = setOfRegionModels[reg];
581 G4int nmod = regModels->NumberOfModels();
582
583 // Calculate energy losses vector
584 std::size_t totBinsLoss = aVector->GetVectorLength();
585 G4double del = 0.0;
586 G4int k0 = 0;
587
588 for(std::size_t j=0; j<totBinsLoss; ++j) {
589 G4double e = aVector->Energy(j);
590
591 // Choose a model of energy losses
592 G4int k = 0;
593 if (nmod > 1) {
594 k = nmod;
595 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
596 do {--k;} while (k>0 && e <= regModels->LowEdgeEnergy(k));
597 //G4cout << "k= " << k << G4endl;
598 if(k > 0 && k != k0) {
599 k0 = k;
600 G4double elow = regModels->LowEdgeEnergy(k);
601 G4double dedx1 =
602 models[regModels->ModelIndex(k-1)]->ComputeDEDX(couple, particle, elow, cut);
603 G4double dedx2 =
604 models[regModels->ModelIndex(k)]->ComputeDEDX(couple, particle, elow, cut);
605 del = (dedx2 > 0.0) ? (dedx1/dedx2 - 1.0)*elow : 0.0;
606 //G4cout << "elow= " << elow
607 // << " dedx1= " << dedx1 << " dedx2= " << dedx2 << G4endl;
608 }
609 }
610 G4double dedx = (1.0 + del/e)*
611 models[regModels->ModelIndex(k)]->ComputeDEDX(couple, particle, e, cut);
612
613 if(2 < verboseLevel) {
614 G4cout << "Material= " << couple->GetMaterial()->GetName()
615 << " E(MeV)= " << e/MeV
616 << " dEdx(MeV/mm)= " << dedx*mm/MeV
617 << " del= " << del*mm/MeV<< " k= " << k
618 << " modelIdx= " << regModels->ModelIndex(k)
619 << G4endl;
620 }
621 dedx = std::max(dedx, 0.0);
622 aVector->PutValue(j, dedx);
623 }
624}
625
626//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627
629 const G4MaterialCutsCouple* couple,
630 G4bool startFromNull,
631 G4EmTableType tType)
632{
633 std::size_t i = couple->GetIndex();
634 G4double cut = (*theCuts)[i];
635 G4double tmax = DBL_MAX;
636
637 G4int reg = 0;
638 if(nRegions > 1 && nEmModels > 1) { reg = idxOfRegionModels[i]; }
639 const G4RegionModels* regModels = setOfRegionModels[reg];
640 G4int nmod = regModels->NumberOfModels();
641 if(1 < verboseLevel) {
642 G4cout << "G4EmModelManager::FillLambdaVector() for "
643 << particle->GetParticleName()
644 << " in " << couple->GetMaterial()->GetName()
645 << " Emin(MeV)= " << aVector->Energy(0)
646 << " Emax(MeV)= " << aVector->GetMaxEnergy()
647 << " cut= " << cut
648 << " Type " << tType
649 << " nmod= " << nmod
650 << G4endl;
651 }
652
653 // Calculate lambda vector
654 std::size_t totBinsLambda = aVector->GetVectorLength();
655 G4double del = 0.0;
656 G4int k0 = 0;
657 G4int k = 0;
658 G4VEmModel* mod = models[regModels->ModelIndex(0)];
659 for(std::size_t j=0; j<totBinsLambda; ++j) {
660
661 G4double e = aVector->Energy(j);
662
663 // Choose a model
664 if (nmod > 1) {
665 k = nmod;
666 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
667 do {--k;} while (k>0 && e <= regModels->LowEdgeEnergy(k));
668 if(k > 0 && k != k0) {
669 k0 = k;
670 G4double elow = regModels->LowEdgeEnergy(k);
671 G4VEmModel* mod1 = models[regModels->ModelIndex(k-1)];
672 G4double xs1 = mod1->CrossSection(couple,particle,elow,cut,tmax);
673 mod = models[regModels->ModelIndex(k)];
674 G4double xs2 = mod->CrossSection(couple,particle,elow,cut,tmax);
675 del = (xs2 > 0.0) ? (xs1/xs2 - 1.0)*elow : 0.0;
676 //G4cout << "New model k=" << k << " E(MeV)= " << e/MeV
677 // << " Elow(MeV)= " << elow/MeV << " del= " << del << G4endl;
678 }
679 }
680 G4double cross = (1.0 + del/e)*mod->CrossSection(couple,particle,e,cut,tmax);
681 if(fIsCrossSectionPrim == tType) { cross *= e; }
682
683 if(j==0 && startFromNull) { cross = 0.0; }
684
685 if(2 < verboseLevel) {
686 G4cout << "FillLambdaVector: " << j << ". e(MeV)= " << e/MeV
687 << " cross(1/mm)= " << cross*mm
688 << " del= " << del*mm << " k= " << k
689 << " modelIdx= " << regModels->ModelIndex(k)
690 << G4endl;
691 }
692 cross = std::max(cross, 0.0);
693 aVector->PutValue(j, cross);
694 }
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
699void G4EmModelManager::DumpModelList(std::ostream& out, G4int verb)
700{
701 if(verb == 0) { return; }
702 for(G4int i=0; i<nRegions; ++i) {
703 G4RegionModels* r = setOfRegionModels[i];
704 const G4Region* reg = r->Region();
705 G4int n = r->NumberOfModels();
706 if(n > 0) {
707 out << " ===== EM models for the G4Region " << reg->GetName()
708 << " ======" << G4endl;
709 for(G4int j=0; j<n; ++j) {
710 G4VEmModel* model = models[r->ModelIndex(j)];
711 G4double emin =
712 std::max(r->LowEdgeEnergy(j),model->LowEnergyActivationLimit());
713 G4double emax =
714 std::min(r->LowEdgeEnergy(j+1),model->HighEnergyActivationLimit());
715 if(emax > emin) {
716 out << std::setw(20);
717 out << model->GetName() << " : Emin="
718 << std::setw(5) << G4BestUnit(emin,"Energy")
719 << " Emax="
720 << std::setw(5) << G4BestUnit(emax,"Energy");
721 G4PhysicsTable* table = model->GetCrossSectionTable();
722 if(table) {
723 std::size_t kk = table->size();
724 for(std::size_t k=0; k<kk; ++k) {
725 const G4PhysicsVector* v = (*table)[k];
726 if(v) {
727 G4int nn = G4int(v->GetVectorLength() - 1);
728 out << " Nbins=" << nn << " "
729 << std::setw(3) << G4BestUnit(v->Energy(0),"Energy")
730 << " - "
731 << std::setw(3) << G4BestUnit(v->Energy(nn),"Energy");
732 break;
733 }
734 }
735 }
737 if(an) { out << " " << an->GetName(); }
738 if(fluoFlag && model->DeexcitationFlag()) {
739 out << " Fluo";
740 }
741 out << G4endl;
742 auto msc = dynamic_cast<G4VMscModel*>(model);
743 if(msc != nullptr) msc->DumpParameters(out);
744 }
745 }
746 }
747 if(1 == nEmModels) { break; }
748 }
749 if(theCutsNew) {
750 out << " ===== Limit on energy threshold has been applied " << G4endl;
751 }
752}
753
754//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4EmTableType
@ fTotal
@ fIsCrossSectionPrim
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4BestUnit(a, b)
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
static G4Electron * Electron()
Definition G4Electron.cc:91
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4int verb)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
G4VEmModel * GetRegionModel(G4int idx, std::size_t index_couple)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
G4int NumberOfRegionModels(std::size_t index_couple) const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
const G4String & GetParticleName() const
void PutValue(const std::size_t index, const G4double value)
G4double GetMaxEnergy() const
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
static G4Positron * Positron()
Definition G4Positron.cc:90
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
const G4String & GetName() const
G4VEmAngularDistribution * GetAngularDistribution()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
virtual void DefineForRegion(const G4Region *)
G4double HighEnergyActivationLimit() const
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4bool DeexcitationFlag() const
const G4String & GetName() const
G4PhysicsTable * GetCrossSectionTable()
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4double LowEnergyActivationLimit() const
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
void DumpParameters(std::ostream &out) const
#define DBL_MAX
Definition templates.hh:62