source: trunk/source/processes/electromagnetic/utils/src/G4EmModelManager.cc@ 1258

Last change on this file since 1258 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 20.7 KB
Line 
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// $Id: G4EmModelManager.cc,v 1.58 2009/10/29 18:07:08 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4EmModelManager
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 07.05.2002
39//
40// Modifications:
41//
42// 23-12-02 V.Ivanchenko change interface in order to move
43// to cut per region
44// 20-01-03 Migrade to cut per region (V.Ivanchenko)
45// 24-01-03 Make models region aware (V.Ivanchenko)
46// 13-02-03 The set of models is defined for region (V.Ivanchenko)
47// 06-03-03 Fix in energy intervals for models (V.Ivanchenko)
48// 13-04-03 Add startFromNull (V.Ivanchenko)
49// 13-05-03 Add calculation of precise range (V.Ivanchenko)
50// 16-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection
51// calculation (V.Ivanchenko)
52// 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
53// 03-11-03 Substitute STL vector for G4RegionModels (V.Ivanchenko)
54// 26-01-04 Fix in energy range conditions (V.Ivanchenko)
55// 24-03-05 Remove check or IsInCharge (V.Ivanchenko)
56// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
57// 18-08-05 Fix cut for e+e- pair production (V.Ivanchenko)
58// 29-11-05 Add protection for arithmetic operations with cut=DBL_MAX (V.Ivanchenko)
59// 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
60// 13-05-06 Add GetModel by index method (VI)
61// 15-03-07 Add maxCutInRange (V.Ivanchenko)
62// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
63// 08-04-08 Fixed and simplified initialisation of G4RegionModel (VI)
64// 03-08-09 Create internal vectors only it is needed (VI)
65//
66// Class Description:
67//
68// It is the unified energy loss process it calculates the continuous
69// energy loss for charged particles using a set of Energy Loss
70// models valid for different energy regions. There are a possibility
71// to create and access to dE/dx and range tables, or to calculate
72// that information on fly.
73// -------------------------------------------------------------------
74//
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78#include "G4EmModelManager.hh"
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
82G4RegionModels::G4RegionModels(G4int nMod, std::vector<G4int>& indx,
83 G4DataVector& lowE, const G4Region* reg)
84{
85 nModelsForRegion = nMod;
86 theListOfModelIndexes = new G4int [nModelsForRegion];
87 lowKineticEnergy = new G4double [nModelsForRegion+1];
88 for (G4int i=0; i<nModelsForRegion; ++i) {
89 theListOfModelIndexes[i] = indx[i];
90 lowKineticEnergy[i] = lowE[i];
91 }
92 lowKineticEnergy[nModelsForRegion] = lowE[nModelsForRegion];
93 theRegion = reg;
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97
98G4RegionModels::~G4RegionModels()
99{
100 delete [] theListOfModelIndexes;
101 delete [] lowKineticEnergy;
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
106#include "G4Step.hh"
107#include "G4ParticleDefinition.hh"
108#include "G4PhysicsVector.hh"
109#include "G4Gamma.hh"
110#include "G4Positron.hh"
111#include "G4MaterialCutsCouple.hh"
112#include "G4ProductionCutsTable.hh"
113#include "G4RegionStore.hh"
114#include "G4Gamma.hh"
115#include "G4Positron.hh"
116#include "G4UnitsTable.hh"
117#include "G4DataVector.hh"
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
121G4EmModelManager::G4EmModelManager():
122 nEmModels(0),
123 nRegions(0),
124 particle(0),
125 verboseLevel(0)
126{
127 maxSubCutInRange = 0.7*mm;
128 models.reserve(4);
129 flucModels.reserve(4);
130 regions.reserve(4);
131 orderOfModels.reserve(4);
132 isUsed.reserve(4);
133 severalModels = true;
134 currRegionModel = 0;
135 currModel = 0;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139
140G4EmModelManager::~G4EmModelManager()
141{
142 verboseLevel = 0; // no verbosity at destruction
143 Clear();
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
148void G4EmModelManager::Clear()
149{
150 if(1 < verboseLevel) {
151 G4cout << "G4EmModelManager::Clear()" << G4endl;
152 }
153 size_t n = setOfRegionModels.size();
154 if(n > 0) {
155 for(size_t i=0; i<n; ++i) {
156 delete setOfRegionModels[i];
157 setOfRegionModels[i] = 0;
158 }
159 }
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163
164void G4EmModelManager::AddEmModel(G4int num, G4VEmModel* p,
165 G4VEmFluctuationModel* fm, const G4Region* r)
166{
167 if(!p) {
168 G4cout << "G4EmModelManager::AddEmModel WARNING: no model defined."
169 << G4endl;
170 return;
171 }
172 models.push_back(p);
173 flucModels.push_back(fm);
174 regions.push_back(r);
175 orderOfModels.push_back(num);
176 isUsed.push_back(0);
177 p->DefineForRegion(r);
178 nEmModels++;
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182
183void G4EmModelManager::UpdateEmModel(const G4String& nam,
184 G4double emin, G4double emax)
185{
186 if (nEmModels > 0) {
187 for(G4int i=0; i<nEmModels; ++i) {
188 if(nam == models[i]->GetName()) {
189 models[i]->SetLowEnergyLimit(emin);
190 models[i]->SetHighEnergyLimit(emax);
191 break;
192 }
193 }
194 }
195 G4cout << "G4EmModelManager::UpdateEmModel WARNING: no model <"
196 << nam << "> is found out"
197 << G4endl;
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
201
202G4VEmModel* G4EmModelManager::GetModel(G4int i, G4bool ver)
203{
204 G4VEmModel* m = 0;
205 if(i >= 0 && i < nEmModels) {m = models[i];}
206 else if(verboseLevel > 0 && ver) {
207 G4cout << "G4EmModelManager::GetModel WARNING: "
208 << "index " << i << " is wrong Nmodels= "
209 << nEmModels;
210 if(particle) G4cout << " for " << particle->GetParticleName();
211 G4cout<< G4endl;
212 }
213 return m;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217
218const G4DataVector*
219G4EmModelManager::Initialise(const G4ParticleDefinition* p,
220 const G4ParticleDefinition* secondaryParticle,
221 G4double minSubRange,
222 G4int val)
223{
224 verboseLevel = val;
225 G4String partname = p->GetParticleName();
226 if(1 < verboseLevel) {
227 G4cout << "G4EmModelManager::Initialise() for "
228 << partname << G4endl;
229 }
230 // Are models defined?
231 if(nEmModels < 1) {
232 G4Exception("G4EmModelManager::Initialise: no model defined for " + partname);
233 }
234
235 particle = p;
236 Clear(); // needed if run is not first
237
238 G4RegionStore* regionStore = G4RegionStore::GetInstance();
239 const G4Region* world =
240 regionStore->GetRegion("DefaultRegionForTheWorld", false);
241
242 // Identify the list of regions with different set of models
243 nRegions = 1;
244 std::vector<const G4Region*> setr;
245 setr.push_back(world);
246 G4bool isWorld = false;
247
248 for (G4int ii=0; ii<nEmModels; ++ii) {
249 const G4Region* r = regions[ii];
250 if ( r == 0 || r == world) {
251 isWorld = true;
252 regions[ii] = world;
253 } else {
254 G4bool newRegion = true;
255 if (nRegions>1) {
256 for (G4int j=1; j<nRegions; ++j) {
257 if ( r == setr[j] ) newRegion = false;
258 }
259 }
260 if (newRegion) {
261 setr.push_back(r);
262 nRegions++;
263 }
264 }
265 }
266 // Are models defined?
267 if(!isWorld) {
268 G4Exception("G4EmModelManager::Initialise: no models defined for " +
269 partname + " in the World volume");
270 }
271
272 G4ProductionCutsTable* theCoupleTable=
273 G4ProductionCutsTable::GetProductionCutsTable();
274 size_t numOfCouples = theCoupleTable->GetTableSize();
275
276 // prepare vectors, shortcut for the case of only 1 model
277 if(nRegions > 1 && nEmModels > 1) {
278 if(numOfCouples > idxOfRegionModels.size()) idxOfRegionModels.resize(numOfCouples);
279 }
280 size_t nr = 1;
281 if(nEmModels > 1) nr = nRegions;
282 if(nr > setOfRegionModels.size()) setOfRegionModels.resize(nr);
283
284 std::vector<G4int> modelAtRegion(nEmModels);
285 std::vector<G4int> modelOrd(nEmModels);
286 G4DataVector eLow(nEmModels+1);
287 G4DataVector eHigh(nEmModels);
288
289 // Order models for regions
290 for (G4int reg=0; reg<nRegions; ++reg) {
291 const G4Region* region = setr[reg];
292 G4int n = 0;
293
294 for (G4int ii=0; ii<nEmModels; ++ii) {
295
296 G4VEmModel* model = models[ii];
297 if ( region == regions[ii] ) {
298
299 G4double tmin = model->LowEnergyLimit();
300 G4double tmax = model->HighEnergyLimit();
301 G4int ord = orderOfModels[ii];
302 G4bool push = true;
303 G4bool insert = false;
304 G4int idx = n;
305
306 if(1 < verboseLevel) {
307 G4cout << "Model #" << ii
308 << " <" << model->GetName() << "> for region <";
309 if (region) G4cout << region->GetName();
310 G4cout << "> "
311 << " tmin(MeV)= " << tmin/MeV
312 << "; tmax(MeV)= " << tmax/MeV
313 << "; order= " << ord
314 << G4endl;
315 }
316
317 if(n > 0) {
318
319 // extend energy range to previous models
320 tmin = std::min(tmin, eHigh[n-1]);
321 tmax = std::max(tmax, eLow[0]);
322 //G4cout << "tmin= " << tmin << " tmax= "
323 // << tmax << " ord= " << ord <<G4endl;
324 // empty energy range
325 if( tmax - tmin <= eV) push = false;
326 // low-energy model
327 else if (tmax == eLow[0]) {
328 push = false;
329 insert = true;
330 idx = 0;
331 // resolve intersections
332 } else if(tmin < eHigh[n-1]) {
333 // compare order
334 for(G4int k=0; k<n; ++k) {
335 // new model has lower application
336 if(ord >= modelOrd[k]) {
337 if(tmin < eHigh[k] && tmin >= eLow[k]) tmin = eHigh[k];
338 if(tmax <= eHigh[k] && tmax > eLow[k]) tmax = eLow[k];
339 if(tmax > eHigh[k] && tmin < eLow[k]) {
340 if(tmax - eHigh[k] > eLow[k] - tmin) tmin = eHigh[k];
341 else tmax = eLow[k];
342 }
343 if( tmax - tmin <= eV) {
344 push = false;
345 break;
346 }
347 }
348 }
349 //G4cout << "tmin= " << tmin << " tmax= "
350 // << tmax << " push= " << push << " idx= " << idx <<G4endl;
351 if(push) {
352 if (tmax == eLow[0]) {
353 push = false;
354 insert = true;
355 idx = 0;
356 // continue resolve intersections
357 } else if(tmin < eHigh[n-1]) {
358 // last energy interval
359 if(tmin > eLow[n-1] && tmax >= eHigh[n-1]) {
360 eHigh[n-1] = tmin;
361 // first energy interval
362 } else if(tmin <= eLow[0] && tmax < eHigh[0]) {
363 eLow[0] = tmax;
364 push = false;
365 insert = true;
366 idx = 0;
367 } else {
368 // find energy interval to replace
369 for(G4int k=0; k<n; ++k) {
370 if(tmin <= eLow[k] && tmax >= eHigh[k]) {
371 push = false;
372 modelAtRegion[k] = ii;
373 modelOrd[k] = ord;
374 isUsed[ii] = 1;
375 }
376 }
377 }
378 }
379 }
380 }
381 }
382 if(insert) {
383 for(G4int k=n-1; k>=idx; --k) {
384 modelAtRegion[k+1] = modelAtRegion[k];
385 modelOrd[k+1] = modelOrd[k];
386 eLow[k+1] = eLow[k];
387 eHigh[k+1] = eHigh[k];
388 }
389 }
390 //G4cout << "push= " << push << " insert= " << insert
391 //<< " idx= " << idx <<G4endl;
392 if (push || insert) {
393 ++n;
394 modelAtRegion[idx] = ii;
395 modelOrd[idx] = ord;
396 eLow[idx] = tmin;
397 eHigh[idx] = tmax;
398 isUsed[ii] = 1;
399 }
400 }
401 }
402 eLow[0] = 0.0;
403 eLow[n] = eHigh[n-1];
404
405 if(1 < verboseLevel) {
406 G4cout << "New G4RegionModels set with " << n << " models for region <";
407 if (region) G4cout << region->GetName();
408 G4cout << "> Elow(MeV)= ";
409 for(G4int ii=0; ii<=n; ++ii) {G4cout << eLow[ii]/MeV << " ";}
410 G4cout << G4endl;
411 }
412 G4RegionModels* rm = new G4RegionModels(n, modelAtRegion, eLow, region);
413 setOfRegionModels[reg] = rm;
414 if(1 == nEmModels) break;
415 }
416
417 currRegionModel = setOfRegionModels[0];
418
419 // Access to materials and build cuts
420 size_t idx = 1;
421 if(secondaryParticle) {
422 if( secondaryParticle == G4Gamma::Gamma() ) idx = 0;
423 else if( secondaryParticle == G4Positron::Positron()) idx = 2;
424 }
425
426 if(numOfCouples > theCuts.size()) {theCuts.resize(numOfCouples);}
427 if(minSubRange < 1.0 && numOfCouples > theSubCuts.size()) {
428 theSubCuts.resize(numOfCouples);
429 }
430 for(size_t i=0; i<numOfCouples; ++i) {
431
432 const G4MaterialCutsCouple* couple =
433 theCoupleTable->GetMaterialCutsCouple(i);
434 const G4Material* material = couple->GetMaterial();
435 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
436
437 G4int reg = 0;
438 if(nRegions > 1 && nEmModels > 1) {
439 reg = nRegions;
440 do {--reg;} while (reg>0 && pcuts != (setr[reg]->GetProductionCuts()));
441 idxOfRegionModels[i] = reg;
442 }
443 if(1 < verboseLevel) {
444 G4cout << "G4EmModelManager::Initialise() for "
445 << material->GetName()
446 << " indexOfCouple= " << i
447 << " indexOfRegion= " << reg
448 << G4endl;
449 }
450
451 G4double cut = (*theCoupleTable->GetEnergyCutsVector(idx))[i];
452 G4double subcut = DBL_MAX;
453 if(secondaryParticle) {
454
455 // compute subcut
456 if( cut < DBL_MAX && minSubRange < 1.0) {
457 subcut = minSubRange*cut;
458 G4double rcut = std::min(minSubRange*pcuts->GetProductionCut(idx),
459 maxSubCutInRange);
460 G4double tcutmax =
461 theCoupleTable->ConvertRangeToEnergy(secondaryParticle,material,rcut);
462 if(tcutmax < subcut) subcut = tcutmax;
463 }
464 }
465
466 G4int nm = setOfRegionModels[reg]->NumberOfModels();
467 for(G4int j=0; j<nm; ++j) {
468
469 G4VEmModel* model = models[setOfRegionModels[reg]->ModelIndex(j)];
470
471 G4double tcutmin = model->MinEnergyCut(particle, couple);
472
473 if(cut < tcutmin) cut = tcutmin;
474 if(subcut < tcutmin) subcut = tcutmin;
475 if(1 < verboseLevel) {
476 G4cout << "The model # " << j
477 << "; tcutmin(MeV)= " << tcutmin/MeV
478 << "; tcut(MeV)= " << cut/MeV
479 << "; tsubcut(MeV)= " << subcut/MeV
480 << " for " << particle->GetParticleName()
481 << G4endl;
482 }
483 }
484 theCuts[i] = cut;
485 if(minSubRange < 1.0) theSubCuts[i] = subcut;
486 }
487
488 // initialize models
489 G4int nn = 0;
490 severalModels = true;
491 for(G4int jj=0; jj<nEmModels; ++jj) {
492 if(1 == isUsed[jj]) {
493 ++nn;
494 currModel = models[jj];
495 currModel->Initialise(particle, theCuts);
496 if(flucModels[jj]) flucModels[jj]->InitialiseMe(particle);
497 }
498 }
499 if(1 == nn) severalModels = false;
500
501 if(1 < verboseLevel) {
502 G4cout << "G4EmModelManager for " << particle->GetParticleName()
503 << " is initialised; nRegions= " << nRegions
504 << G4endl;
505 }
506
507 return &theCuts;
508}
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
511
512void G4EmModelManager::FillDEDXVector(G4PhysicsVector* aVector,
513 const G4MaterialCutsCouple* couple,
514 G4EmTableType tType)
515{
516 size_t i = couple->GetIndex();
517 G4double cut = theCuts[i];
518 G4double emin = 0.0;
519
520 if(fTotal == tType) cut = DBL_MAX;
521 else if(fSubRestricted == tType) {
522 emin = cut;
523 if(theSubCuts.size() > 0) emin = theSubCuts[i];
524 }
525
526 if(1 < verboseLevel) {
527 G4cout << "G4EmModelManager::FillDEDXVector() for "
528 << couple->GetMaterial()->GetName()
529 << " cut(MeV)= " << cut
530 << " emin(MeV)= " << emin
531 << " Type " << tType
532 << " for " << particle->GetParticleName()
533 << G4endl;
534 }
535
536 G4int reg = 0;
537 if(nRegions > 1 && nEmModels > 1) reg = idxOfRegionModels[i];
538 const G4RegionModels* regModels = setOfRegionModels[reg];
539 G4int nmod = regModels->NumberOfModels();
540
541 // Calculate energy losses vector
542
543 //G4cout << "nmod= " << nmod << G4endl;
544 size_t totBinsLoss = aVector->GetVectorLength();
545 for(size_t j=0; j<totBinsLoss; ++j) {
546
547 G4double e = aVector->Energy(j);
548 G4double del = 0.0;
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) {
557 G4double elow = regModels->LowEdgeEnergy(k);
558 G4double dedx1 = ComputeDEDX(models[regModels->ModelIndex(k-1)],
559 couple,elow,cut,emin);
560 G4double dedx2 = ComputeDEDX(models[regModels->ModelIndex(k)],
561 couple,elow,cut,emin);
562 del = (dedx1 - dedx2)*elow/e;
563 //G4cout << "elow= " << elow
564 // << " dedx1= " << dedx1 << " dedx2= " << dedx2 << G4endl;
565 }
566 }
567 G4double dedx = ComputeDEDX(models[regModels->ModelIndex(k)],
568 couple,e,cut,emin) + del;
569
570 if(dedx < 0.0) dedx = 0.0;
571 if(2 < verboseLevel) {
572 G4cout << "Material= " << couple->GetMaterial()->GetName()
573 << " E(MeV)= " << e/MeV
574 << " dEdx(MeV/mm)= " << dedx*mm/MeV
575 << " del= " << del*mm/MeV<< " k= " << k
576 << " modelIdx= " << regModels->ModelIndex(k)
577 << G4endl;
578 }
579 aVector->PutValue(j, dedx);
580 }
581}
582
583//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584
585void G4EmModelManager::FillLambdaVector(G4PhysicsVector* aVector,
586 const G4MaterialCutsCouple* couple,
587 G4bool startFromNull,
588 G4EmTableType tType)
589{
590 size_t i = couple->GetIndex();
591 G4double cut = theCuts[i];
592 G4double tmax = DBL_MAX;
593 if (fSubRestricted == tType) {
594 tmax = cut;
595 if(theSubCuts.size() > 0) cut = theSubCuts[i];
596 }
597
598 G4int reg = 0;
599 if(nRegions > 1 && nEmModels > 1) reg = idxOfRegionModels[i];
600 const G4RegionModels* regModels = setOfRegionModels[reg];
601 G4int nmod = regModels->NumberOfModels();
602 if(1 < verboseLevel) {
603 G4cout << "G4EmModelManager::FillLambdaVector() for "
604 << particle->GetParticleName()
605 << " in " << couple->GetMaterial()->GetName()
606 << " Ecut(MeV)= " << cut
607 << " Emax(MeV)= " << tmax
608 << " Type " << tType
609 << " nmod= " << nmod
610 << G4endl;
611 }
612
613
614 // Calculate lambda vector
615 size_t totBinsLambda = aVector->GetVectorLength();
616 for(size_t j=0; j<totBinsLambda; ++j) {
617
618 G4double e = aVector->Energy(j);
619
620 G4double del = 0.0;
621
622 // Choose a model
623 G4int k = 0;
624 G4VEmModel* mod = models[regModels->ModelIndex(0)];
625 if (nmod > 1) {
626 k = nmod;
627 do {--k;} while (k>0 && e <= regModels->LowEdgeEnergy(k));
628 if(k > 0) {
629 G4double elow = regModels->LowEdgeEnergy(k);
630 G4VEmModel* m = models[regModels->ModelIndex(k-1)];
631 G4double xs1 = m->CrossSection(couple,particle,elow,cut,tmax);
632 mod = models[regModels->ModelIndex(k)];
633 G4double xs2 = mod->CrossSection(couple,particle,elow,cut,tmax);
634 del = (xs1 - xs2)*elow/e;
635 }
636 }
637 G4double cross = mod->CrossSection(couple,particle,e,cut,tmax) + del;
638
639 if(j==0 && startFromNull) cross = 0.0;
640
641 if(2 < verboseLevel) {
642 G4cout << "FillLambdaVector: " << j << ". e(MeV)= " << e/MeV
643 << " cross(1/mm)= " << cross*mm
644 << " del= " << del*mm << " k= " << k
645 << " modelIdx= " << regModels->ModelIndex(k)
646 << G4endl;
647 }
648 if(cross < 0.0) cross = 0.0;
649
650 aVector->PutValue(j, cross);
651 }
652}
653
654//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
655
656void G4EmModelManager::DumpModelList(G4int verb)
657{
658 if(verb == 0) return;
659 for(G4int i=0; i<nRegions; ++i) {
660 G4RegionModels* r = setOfRegionModels[i];
661 const G4Region* reg = r->Region();
662 G4int n = r->NumberOfModels();
663 if(n > 0) {
664 G4cout << " ===== EM models for the G4Region " << reg->GetName()
665 << " ======" << G4endl;;
666 for(G4int j=0; j<n; ++j) {
667 const G4VEmModel* m = models[r->ModelIndex(j)];
668 G4cout << std::setw(20);
669 G4cout << m->GetName() << " : Emin= "
670 << std::setw(10) << G4BestUnit(r->LowEdgeEnergy(j),"Energy")
671 << " Emax= "
672 << G4BestUnit(r->LowEdgeEnergy(j+1),"Energy")
673 << G4endl;
674 }
675 }
676 if(1 == nEmModels) break;
677 }
678}
679
680//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.