Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4VEmProcess.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4VEmProcess
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 01.10.2003
38//
39// Modifications:
40// 30-06-04 make it to be pure discrete process (V.Ivanchenko)
41// 30-09-08 optimise integral option (V.Ivanchenko)
42// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
43// 11-03-05 Shift verbose level by 1, add applyCuts and killPrimary flags (VI)
44// 14-03-05 Update logic PostStepDoIt (V.Ivanchenko)
45// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
46// 18-04-05 Use G4ParticleChangeForGamma (V.Ivanchenko)
47// 25-07-05 Add protection: integral mode only for charged particles (VI)
48// 04-09-05 default lambdaFactor 0.8 (V.Ivanchenko)
49// 11-01-06 add A to parameters of ComputeCrossSectionPerAtom (VI)
50// 12-09-06 add SetModel() (mma)
51// 12-04-07 remove double call to Clear model manager (V.Ivanchenko)
52// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
53// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
54// 17-02-10 Added pointer currentParticle (VI)
55// 30-05-12 allow Russian roulette, brem splitting (D. Sawkey)
56//
57// Class Description:
58//
59
60// -------------------------------------------------------------------
61//
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
65#include "G4VEmProcess.hh"
67#include "G4SystemOfUnits.hh"
68#include "G4ProcessManager.hh"
69#include "G4LossTableManager.hh"
70#include "G4LossTableBuilder.hh"
71#include "G4Step.hh"
73#include "G4VEmModel.hh"
74#include "G4DataVector.hh"
75#include "G4PhysicsTable.hh"
76#include "G4PhysicsLogVector.hh"
77#include "G4VParticleChange.hh"
79#include "G4Region.hh"
80#include "G4Gamma.hh"
81#include "G4Electron.hh"
82#include "G4Positron.hh"
84#include "G4EmBiasingManager.hh"
85#include "G4GenericIon.hh"
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90 G4VDiscreteProcess(name, type),
91 secondaryParticle(0),
92 buildLambdaTable(true),
93 numberOfModels(0),
94 theLambdaTable(0),
95 theLambdaTablePrim(0),
96 theDensityFactor(0),
97 theDensityIdx(0),
98 integral(false),
99 applyCuts(false),
100 startFromNull(false),
101 splineFlag(true),
102 currentModel(0),
103 particle(0),
104 currentParticle(0),
105 currentCouple(0)
106{
108
109 // Size of tables assuming spline
110 minKinEnergy = 0.1*keV;
111 maxKinEnergy = 10.0*TeV;
112 nLambdaBins = 77;
113 minKinEnergyPrim = DBL_MAX;
114
115 // default lambda factor
116 lambdaFactor = 0.8;
117
118 // default limit on polar angle
119 polarAngleLimit = 0.0;
120 biasFactor = 1.0;
121
122 // particle types
123 theGamma = G4Gamma::Gamma();
124 theElectron = G4Electron::Electron();
125 thePositron = G4Positron::Positron();
126
128 secParticles.reserve(5);
129
130 preStepLambda = 0.0;
131 mfpKinEnergy = DBL_MAX;
132
133 modelManager = new G4EmModelManager();
134 biasManager = 0;
135 biasFlag = false;
136 weightFlag = false;
137 (G4LossTableManager::Instance())->Register(this);
138 warn = 0;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144{
145 if(1 < verboseLevel) {
146 G4cout << "G4VEmProcess destruct " << GetProcessName()
147 << " " << this << " " << theLambdaTable <<G4endl;
148 }
149 Clear();
150 if(theLambdaTable) {
151 theLambdaTable->clearAndDestroy();
152 delete theLambdaTable;
153 }
154 if(theLambdaTablePrim) {
155 theLambdaTablePrim->clearAndDestroy();
156 delete theLambdaTablePrim;
157 }
158 delete modelManager;
159 delete biasManager;
160 (G4LossTableManager::Instance())->DeRegister(this);
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
165void G4VEmProcess::Clear()
166{
167 currentCouple = 0;
168 preStepLambda = 0.0;
169 mfpKinEnergy = DBL_MAX;
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173
175 const G4Material*)
176{
177 return 0.0;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
183 const G4Region* region)
184{
185 G4VEmFluctuationModel* fm = 0;
186 modelManager->AddEmModel(order, p, fm, region);
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
191
193{
194 ++warn;
195 if(warn < 10) {
196 G4cout << "### G4VEmProcess::SetModel is obsolete method and will be "
197 << "removed for the next release." << G4endl;
198 G4cout << " Please, use SetEmModel" << G4endl;
199 }
200 G4int n = emModels.size();
201 if(index >= n) { for(G4int i=n; i<=index; ++i) {emModels.push_back(0);} }
202 emModels[index] = p;
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
208{
209 if(warn < 10) {
210 G4cout << "### G4VEmProcess::Model is obsolete method and will be "
211 << "removed for the next release." << G4endl;
212 G4cout << " Please, use EmModel" << G4endl;
213 }
214 G4VEmModel* p = 0;
215 if(index >= 0 && index < G4int(emModels.size())) { p = emModels[index]; }
216 return p;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220
222{
223 G4int n = emModels.size();
224 if(index >= n) { for(G4int i=n; i<=index; ++i) {emModels.push_back(0);} }
225 emModels[index] = p;
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
231{
232 G4VEmModel* p = 0;
233 if(index >= 0 && index < G4int(emModels.size())) { p = emModels[index]; }
234 return p;
235}
236
237//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238
240 G4double emin, G4double emax)
241{
242 modelManager->UpdateEmModel(nam, emin, emax);
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246
248{
249 return modelManager->GetModel(idx, ver);
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253
255{
256 if(!particle) { SetParticle(&part); }
257
258 if(part.GetParticleType() == "nucleus" &&
259 part.GetParticleSubType() == "generic") {
260
261 G4String pname = part.GetParticleName();
262 if(pname != "deuteron" && pname != "triton" &&
263 pname != "alpha" && pname != "He3" &&
264 pname != "alpha+" && pname != "helium" &&
265 pname != "hydrogen") {
266
267 particle = G4GenericIon::GenericIon();
268 }
269 }
270
271 if(1 < verboseLevel) {
272 G4cout << "G4VEmProcess::PreparePhysicsTable() for "
273 << GetProcessName()
274 << " and particle " << part.GetParticleName()
275 << " local particle " << particle->GetParticleName()
276 << G4endl;
277 }
278
281
282 man->PreparePhysicsTable(&part, this);
283
284 if(particle == &part) {
285 Clear();
286 InitialiseProcess(particle);
287
288 const G4ProductionCutsTable* theCoupleTable=
290 size_t n = theCoupleTable->GetTableSize();
291
292 theEnergyOfCrossSectionMax.resize(n, 0.0);
293 theCrossSectionMax.resize(n, DBL_MAX);
294
295 // initialisation of models
296 numberOfModels = modelManager->NumberOfModels();
297 for(G4int i=0; i<numberOfModels; ++i) {
298 G4VEmModel* mod = modelManager->GetModel(i);
299 if(0 == i) { currentModel = mod; }
300 mod->SetPolarAngleLimit(polarAngleLimit);
301 if(mod->HighEnergyLimit() > maxKinEnergy) {
302 mod->SetHighEnergyLimit(maxKinEnergy);
303 }
304 }
305
306 if(man->AtomDeexcitation()) { modelManager->SetFluoFlag(true); }
307 theCuts = modelManager->Initialise(particle,secondaryParticle,
308 2.,verboseLevel);
309 theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
310 theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
311 theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
312
313 // prepare tables
314 if(buildLambdaTable){
315 theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
316 bld->InitialiseBaseMaterials(theLambdaTable);
317 }
318 // high energy table
319 if(minKinEnergyPrim < maxKinEnergy){
320 theLambdaTablePrim =
322 bld->InitialiseBaseMaterials(theLambdaTablePrim);
323 }
324 // forced biasing
325 if(biasManager) {
326 biasManager->Initialise(part,GetProcessName(),verboseLevel);
327 biasFlag = false;
328 }
329 }
330 theDensityFactor = bld->GetDensityFactors();
331 theDensityIdx = bld->GetCoupleIndexes();
332}
333
334//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
335
337{
338 G4String num = part.GetParticleName();
339 if(1 < verboseLevel) {
340 G4cout << "G4VEmProcess::BuildPhysicsTable() for "
341 << GetProcessName()
342 << " and particle " << num
343 << " buildLambdaTable= " << buildLambdaTable
344 << G4endl;
345 }
346
348
349 if(buildLambdaTable || minKinEnergyPrim < maxKinEnergy) {
350 BuildLambdaTable();
351 }
352
353 // explicitly defined printout by particle name
354 if(1 < verboseLevel ||
355 (0 < verboseLevel && (num == "gamma" || num == "e-" ||
356 num == "e+" || num == "mu+" ||
357 num == "mu-" || num == "proton"||
358 num == "pi+" || num == "pi-" ||
359 num == "kaon+" || num == "kaon-" ||
360 num == "alpha" || num == "anti_proton" ||
361 num == "GenericIon")))
362 {
363 particle = &part;
365 }
366
367 if(1 < verboseLevel) {
368 G4cout << "G4VEmProcess::BuildPhysicsTable() done for "
369 << GetProcessName()
370 << " and particle " << num
371 << G4endl;
372 }
373}
374
375//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
376
377void G4VEmProcess::BuildLambdaTable()
378{
379 if(1 < verboseLevel) {
380 G4cout << "G4EmProcess::BuildLambdaTable() for process "
381 << GetProcessName() << " and particle "
382 << particle->GetParticleName() << " " << this
383 << G4endl;
384 }
385
386 // Access to materials
387 const G4ProductionCutsTable* theCoupleTable=
389 size_t numOfCouples = theCoupleTable->GetTableSize();
390
391 G4LossTableBuilder* bld = (G4LossTableManager::Instance())->GetTableBuilder();
392
393 G4PhysicsLogVector* aVector = 0;
394 G4PhysicsLogVector* bVector = 0;
395 G4PhysicsLogVector* aVectorPrim = 0;
396 G4PhysicsLogVector* bVectorPrim = 0;
397
398 G4double scale = 1.0;
399 G4double emax1 = maxKinEnergy;
400 if(startFromNull || minKinEnergyPrim < maxKinEnergy ) {
401 scale = std::log(maxKinEnergy/minKinEnergy);
402 if(minKinEnergyPrim < maxKinEnergy) { emax1 = minKinEnergyPrim; }
403 }
404
405 for(size_t i=0; i<numOfCouples; ++i) {
406
407 if (bld->GetFlag(i)) {
408
409 // create physics vector and fill it
410 const G4MaterialCutsCouple* couple =
411 theCoupleTable->GetMaterialCutsCouple(i);
412
413 // build main table
414 if(buildLambdaTable) {
415 delete (*theLambdaTable)[i];
416
417 G4bool startNull = startFromNull;
418 // if start from zero then change the scale
419 if(startFromNull || minKinEnergyPrim < maxKinEnergy) {
420 G4double emin = MinPrimaryEnergy(particle,couple->GetMaterial());
421 if(emin < minKinEnergy) {
422 emin = minKinEnergy;
423 startNull = false;
424 }
425 G4double emax = emax1;
426 if(emax <= emin) { emax = 2*emin; }
427 G4int bin =
428 G4lrint(nLambdaBins*std::log(emax/emin)/scale);
429 if(bin < 3) { bin = 3; }
430 aVector = new G4PhysicsLogVector(emin, emax, bin);
431
432 // start not from zero
433 } else if(!bVector) {
434 aVector =
435 new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nLambdaBins);
436 bVector = aVector;
437 } else {
438 aVector = new G4PhysicsLogVector(*bVector);
439 }
440 aVector->SetSpline(splineFlag);
441 modelManager->FillLambdaVector(aVector, couple, startNull);
442 if(splineFlag) { aVector->FillSecondDerivatives(); }
443 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector);
444 }
445 // build high energy table
446 if(minKinEnergyPrim < maxKinEnergy) {
447 delete (*theLambdaTablePrim)[i];
448
449 // start not from zero
450 if(!bVectorPrim) {
451 G4int bin =
452 G4lrint(nLambdaBins*std::log(maxKinEnergy/minKinEnergyPrim)/scale);
453 if(bin < 3) { bin = 3; }
454 aVectorPrim =
455 new G4PhysicsLogVector(minKinEnergyPrim, maxKinEnergy, bin);
456 bVectorPrim = aVectorPrim;
457 } else {
458 aVectorPrim = new G4PhysicsLogVector(*bVectorPrim);
459 }
460 // always use spline
461 aVectorPrim->SetSpline(true);
462 modelManager->FillLambdaVector(aVectorPrim, couple, false,
464 aVectorPrim->FillSecondDerivatives();
465 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTablePrim, i, aVectorPrim);
466 }
467 }
468 }
469
470 if(buildLambdaTable) { FindLambdaMax(); }
471
472 if(1 < verboseLevel) {
473 G4cout << "Lambda table is built for "
474 << particle->GetParticleName()
475 << G4endl;
476 }
477}
478
479//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480
482{
483 if(verboseLevel > 0) {
484 G4cout << G4endl << GetProcessName() << ": for "
485 << particle->GetParticleName();
486 if(integral) { G4cout << ", integral: 1 "; }
487 if(applyCuts) { G4cout << ", applyCuts: 1 "; }
488 G4cout << " SubType= " << GetProcessSubType();;
489 if(biasFactor != 1.0) { G4cout << " BiasingFactor= " << biasFactor; }
490 G4cout << G4endl;
491 if(buildLambdaTable) {
492 size_t length = theLambdaTable->length();
493 for(size_t i=0; i<length; ++i) {
494 G4PhysicsVector* v = (*theLambdaTable)[i];
495 if(v) {
496 G4cout << " Lambda table from "
497 << G4BestUnit(minKinEnergy,"Energy")
498 << " to "
499 << G4BestUnit(v->GetMaxEnergy(),"Energy")
500 << " in " << v->GetVectorLength()-1
501 << " bins, spline: "
502 << splineFlag
503 << G4endl;
504 break;
505 }
506 }
507 }
508 if(minKinEnergyPrim < maxKinEnergy) {
509 size_t length = theLambdaTablePrim->length();
510 for(size_t i=0; i<length; ++i) {
511 G4PhysicsVector* v = (*theLambdaTablePrim)[i];
512 if(v) {
513 G4cout << " LambdaPrime table from "
514 << G4BestUnit(minKinEnergyPrim,"Energy")
515 << " to "
516 << G4BestUnit(maxKinEnergy,"Energy")
517 << " in " << v->GetVectorLength()-1
518 << " bins "
519 << G4endl;
520 break;
521 }
522 }
523 }
524 PrintInfo();
525 modelManager->DumpModelList(verboseLevel);
526 }
527
528 if(verboseLevel > 2 && buildLambdaTable) {
529 G4cout << " LambdaTable address= " << theLambdaTable << G4endl;
530 if(theLambdaTable) { G4cout << (*theLambdaTable) << G4endl; }
531 }
532}
533
534//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
535
537{
538 // reset parameters for the new track
539 currentParticle = track->GetParticleDefinition();
541 //currentInteractionLength = -1.0;
542 // theInitialNumberOfInteractionLength=-1.0;
543 mfpKinEnergy = DBL_MAX;
544
545 // forced biasing only for primary particles
546 if(biasManager) {
547 if(0 == track->GetParentID()) {
548 // primary particle
549 biasFlag = true;
550 biasManager->ResetForcedInteraction();
551 }
552 }
553}
554
555//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
556
558 const G4Track& track,
559 G4double previousStepSize,
561{
563 G4double x = DBL_MAX;
564
565 preStepKinEnergy = track.GetKineticEnergy();
566 DefineMaterial(track.GetMaterialCutsCouple());
567 SelectModel(preStepKinEnergy, currentCoupleIndex);
568
569 if(!currentModel->IsActive(preStepKinEnergy)) { return x; }
570
571 // forced biasing only for primary particles
572 if(biasManager) {
573 if(0 == track.GetParentID()) {
574 if(biasFlag && biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
575 return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
576 }
577 }
578 }
579
580 // compute mean free path
581 if(preStepKinEnergy < mfpKinEnergy) {
582 if (integral) { ComputeIntegralLambda(preStepKinEnergy); }
583 else { preStepLambda = GetCurrentLambda(preStepKinEnergy); }
584
585 // zero cross section
586 if(preStepLambda <= 0.0) {
589 }
590 }
591
592 // non-zero cross section
593 if(preStepLambda > 0.0) {
594
596
597 // beggining of tracking (or just after DoIt of this process)
599
600 } else if(currentInteractionLength < DBL_MAX) {
601
602 // subtract NumberOfInteractionLengthLeft using previous step
604 //SubtractNumberOfInteractionLengthLeft(previousStepSize);
607 //theNumberOfInteractionLengthLeft = perMillion;
608 }
609 }
610
611 // new mean free path and step limit for the next step
612 currentInteractionLength = 1.0/preStepLambda;
614#ifdef G4VERBOSE
615 if (verboseLevel>2){
616 G4cout << "G4VEmProcess::PostStepGetPhysicalInteractionLength ";
617 G4cout << "[ " << GetProcessName() << "]" << G4endl;
618 G4cout << " for " << currentParticle->GetParticleName()
619 << " in Material " << currentMaterial->GetName()
620 << " Ekin(MeV)= " << preStepKinEnergy/MeV
621 <<G4endl;
622 G4cout << " MeanFreePath = " << currentInteractionLength/cm << "[cm]"
623 << " InteractionLength= " << x/cm <<"[cm] " <<G4endl;
624 }
625#endif
626 }
627 return x;
628}
629
630//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
631
633 const G4Step& step)
634{
635 // In all cases clear number of interaction lengths
637 mfpKinEnergy = DBL_MAX;
638
640
641 // Do not make anything if particle is stopped, the annihilation then
642 // should be performed by the AtRestDoIt!
643 if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
644
645 G4double finalT = track.GetKineticEnergy();
646
647 // forced process - should happen only once per track
648 if(biasFlag) {
649 if(biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
650 biasFlag = false;
651 }
652 }
653
654 // Integral approach
655 if (integral) {
656 G4double lx = GetLambda(finalT, currentCouple);
657 if(preStepLambda<lx && 1 < verboseLevel) {
658 G4cout << "WARNING: for " << currentParticle->GetParticleName()
659 << " and " << GetProcessName()
660 << " E(MeV)= " << finalT/MeV
661 << " preLambda= " << preStepLambda << " < "
662 << lx << " (postLambda) "
663 << G4endl;
664 }
665
666 if(preStepLambda*G4UniformRand() > lx) {
668 return &fParticleChange;
669 }
670 }
671
672 SelectModel(finalT, currentCoupleIndex);
673 if(!currentModel->IsActive(finalT)) { return &fParticleChange; }
674
675 // define new weight for primary and secondaries
677 if(weightFlag) {
678 weight /= biasFactor;
680 }
681
682 /*
683 if(0 < verboseLevel) {
684 G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
685 << finalT/MeV
686 << " MeV; model= (" << currentModel->LowEnergyLimit()
687 << ", " << currentModel->HighEnergyLimit() << ")"
688 << G4endl;
689 }
690 */
691
692 // sample secondaries
693 secParticles.clear();
694 currentModel->SampleSecondaries(&secParticles,
695 currentCouple,
696 track.GetDynamicParticle(),
697 (*theCuts)[currentCoupleIndex]);
698
699 // bremsstrahlung splitting or Russian roulette
700 if(biasManager) {
701 if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
702 G4double eloss = 0.0;
703 weight *= biasManager->ApplySecondaryBiasing(secParticles,
704 track, currentModel,
706 eloss, currentCoupleIndex,
707 (*theCuts)[currentCoupleIndex],
708 step.GetPostStepPoint()->GetSafety());
709 if(eloss > 0.0) {
712 }
713 }
714 }
715
716 // save secondaries
717 G4int num = secParticles.size();
718 if(num > 0) {
719
722
723 for (G4int i=0; i<num; ++i) {
724 if (secParticles[i]) {
725 G4DynamicParticle* dp = secParticles[i];
727 G4double e = dp->GetKineticEnergy();
728 G4bool good = true;
729 if(applyCuts) {
730 if (p == theGamma) {
731 if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
732
733 } else if (p == theElectron) {
734 if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
735
736 } else if (p == thePositron) {
737 if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
738 e < (*theCutsPositron)[currentCoupleIndex]) {
739 good = false;
740 e += 2.0*electron_mass_c2;
741 }
742 }
743 // added secondary if it is good
744 }
745 if (good) {
746 G4Track* t = new G4Track(dp, track.GetGlobalTime(), track.GetPosition());
748 t->SetWeight(weight);
750 //G4cout << "Secondary(post step) has weight " << t->GetWeight()
751 // << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
752 } else {
753 delete dp;
754 edep += e;
755 }
756 }
757 }
759 }
760
763 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
766 }
767
768 // ClearNumberOfInteractionLengthLeft();
769 return &fParticleChange;
770}
771
772//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773
775 const G4String& directory,
776 G4bool ascii)
777{
778 G4bool yes = true;
779
780 if ( theLambdaTable && part == particle) {
781 const G4String name =
782 GetPhysicsTableFileName(part,directory,"Lambda",ascii);
783 yes = theLambdaTable->StorePhysicsTable(name,ascii);
784
785 if ( yes ) {
786 G4cout << "Physics table is stored for " << particle->GetParticleName()
787 << " and process " << GetProcessName()
788 << " in the directory <" << directory
789 << "> " << G4endl;
790 } else {
791 G4cout << "Fail to store Physics Table for "
792 << particle->GetParticleName()
793 << " and process " << GetProcessName()
794 << " in the directory <" << directory
795 << "> " << G4endl;
796 }
797 }
798 if ( theLambdaTablePrim && part == particle) {
799 const G4String name =
800 GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
801 yes = theLambdaTablePrim->StorePhysicsTable(name,ascii);
802
803 if ( yes ) {
804 G4cout << "Physics table prim is stored for " << particle->GetParticleName()
805 << " and process " << GetProcessName()
806 << " in the directory <" << directory
807 << "> " << G4endl;
808 } else {
809 G4cout << "Fail to store Physics Table Prim for "
810 << particle->GetParticleName()
811 << " and process " << GetProcessName()
812 << " in the directory <" << directory
813 << "> " << G4endl;
814 }
815 }
816 return yes;
817}
818
819//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
820
822 const G4String& directory,
823 G4bool ascii)
824{
825 if(1 < verboseLevel) {
826 G4cout << "G4VEmProcess::RetrievePhysicsTable() for "
827 << part->GetParticleName() << " and process "
828 << GetProcessName() << G4endl;
829 }
830 G4bool yes = true;
831
832 if((!buildLambdaTable && minKinEnergyPrim > maxKinEnergy)
833 || particle != part) { return yes; }
834
835 const G4String particleName = part->GetParticleName();
836 G4String filename;
837
838 if(buildLambdaTable) {
839 filename = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
841 filename,ascii);
842 if ( yes ) {
843 if (0 < verboseLevel) {
844 G4cout << "Lambda table for " << particleName
845 << " is Retrieved from <"
846 << filename << ">"
847 << G4endl;
848 }
849 if((G4LossTableManager::Instance())->SplineFlag()) {
850 size_t n = theLambdaTable->length();
851 for(size_t i=0; i<n; ++i) {
852 if((* theLambdaTable)[i]) {
853 (* theLambdaTable)[i]->SetSpline(true);
854 }
855 }
856 }
857 } else {
858 if (1 < verboseLevel) {
859 G4cout << "Lambda table for " << particleName << " in file <"
860 << filename << "> is not exist"
861 << G4endl;
862 }
863 }
864 }
865 if(minKinEnergyPrim < maxKinEnergy) {
866 filename = GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
867 yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTablePrim,
868 filename,ascii);
869 if ( yes ) {
870 if (0 < verboseLevel) {
871 G4cout << "Lambda table prim for " << particleName
872 << " is Retrieved from <"
873 << filename << ">"
874 << G4endl;
875 }
876 if((G4LossTableManager::Instance())->SplineFlag()) {
877 size_t n = theLambdaTablePrim->length();
878 for(size_t i=0; i<n; ++i) {
879 if((* theLambdaTablePrim)[i]) {
880 (* theLambdaTablePrim)[i]->SetSpline(true);
881 }
882 }
883 }
884 } else {
885 if (1 < verboseLevel) {
886 G4cout << "Lambda table prim for " << particleName << " in file <"
887 << filename << "> is not exist"
888 << G4endl;
889 }
890 }
891 }
892
893 return yes;
894}
895
896//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
897
900 const G4MaterialCutsCouple* couple)
901{
902 // Cross section per atom is calculated
903 DefineMaterial(couple);
904 G4double cross = 0.0;
905 if(theLambdaTable) {
906 cross = (*theDensityFactor)[currentCoupleIndex]*
907 (((*theLambdaTable)[basedCoupleIndex])->Value(kineticEnergy));
908 } else {
909 SelectModel(kineticEnergy, currentCoupleIndex);
910 cross = currentModel->CrossSectionPerVolume(currentMaterial,
911 currentParticle,kineticEnergy);
912 }
913
914 if(cross < 0.0) { cross = 0.0; }
915 return cross;
916}
917
918//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
919
921 G4double,
923{
925 return G4VEmProcess::MeanFreePath(track);
926}
927
928//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
929
931{
932 DefineMaterial(track.GetMaterialCutsCouple());
933 preStepLambda = GetCurrentLambda(track.GetKineticEnergy());
934 G4double x = DBL_MAX;
935 if(0.0 < preStepLambda) { x = 1.0/preStepLambda; }
936 return x;
937}
938
939//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
940
943 G4double Z, G4double A, G4double cut)
944{
945 SelectModel(kineticEnergy, currentCoupleIndex);
946 G4double x = 0.0;
947 if(currentModel) {
948 x = currentModel->ComputeCrossSectionPerAtom(currentParticle,kineticEnergy,
949 Z,A,cut);
950 }
951 return x;
952}
953
954//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
955
956void G4VEmProcess::FindLambdaMax()
957{
958 if(1 < verboseLevel) {
959 G4cout << "### G4VEmProcess::FindLambdaMax: "
960 << particle->GetParticleName()
961 << " and process " << GetProcessName() << G4endl;
962 }
963 size_t n = theLambdaTable->length();
964 G4PhysicsVector* pv;
965 G4double e, ss, emax, smax;
966
967 size_t i;
968
969 // first loop on existing vectors
970 for (i=0; i<n; ++i) {
971 pv = (*theLambdaTable)[i];
972 if(pv) {
973 size_t nb = pv->GetVectorLength();
974 emax = DBL_MAX;
975 smax = 0.0;
976 if(nb > 0) {
977 for (size_t j=0; j<nb; ++j) {
978 e = pv->Energy(j);
979 ss = (*pv)(j);
980 if(ss > smax) {
981 smax = ss;
982 emax = e;
983 }
984 }
985 }
986 theEnergyOfCrossSectionMax[i] = emax;
987 theCrossSectionMax[i] = smax;
988 if(1 < verboseLevel) {
989 G4cout << "For " << particle->GetParticleName()
990 << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
991 << " lambda= " << smax << G4endl;
992 }
993 }
994 }
995 // second loop using base materials
996 for (i=0; i<n; ++i) {
997 pv = (*theLambdaTable)[i];
998 if(!pv){
999 G4int j = (*theDensityIdx)[i];
1000 theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
1001 theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
1002 }
1003 }
1004}
1005
1006//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1007
1010{
1011 G4PhysicsVector* v =
1012 new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nLambdaBins);
1013 v->SetSpline((G4LossTableManager::Instance())->SplineFlag());
1014 return v;
1015}
1016
1017//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1018
1020{
1021 const G4Element* elm = 0;
1022 if(currentModel) {elm = currentModel->GetCurrentElement(); }
1023 return elm;
1024}
1025
1026//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1027
1029{
1030 if(f > 0.0) {
1031 biasFactor = f;
1032 weightFlag = flag;
1033 if(1 < verboseLevel) {
1034 G4cout << "### SetCrossSectionBiasingFactor: for "
1035 << particle->GetParticleName()
1036 << " and process " << GetProcessName()
1037 << " biasFactor= " << f << " weightFlag= " << flag
1038 << G4endl;
1039 }
1040 }
1041}
1042
1043//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1044
1045void
1047 G4bool flag)
1048{
1049 if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1050 if(1 < verboseLevel) {
1051 G4cout << "### ActivateForcedInteraction: for "
1052 << particle->GetParticleName()
1053 << " and process " << GetProcessName()
1054 << " length(mm)= " << length/mm
1055 << " in G4Region <" << r
1056 << "> weightFlag= " << flag
1057 << G4endl;
1058 }
1059 weightFlag = flag;
1060 biasManager->ActivateForcedInteraction(length, r);
1061}
1062
1063//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1064
1065void
1067 G4double factor,
1068 G4double energyLimit)
1069{
1070 if (0.0 <= factor) {
1071
1072 // Range cut can be applied only for e-
1073 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1074 { return; }
1075
1076 if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1077 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1078 if(1 < verboseLevel) {
1079 G4cout << "### ActivateSecondaryBiasing: for "
1080 << " process " << GetProcessName()
1081 << " factor= " << factor
1082 << " in G4Region <" << region
1083 << "> energyLimit(MeV)= " << energyLimit/MeV
1084 << G4endl;
1085 }
1086 }
1087}
1088
1089//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1090
1092{
1093 nLambdaBins = G4lrint(nLambdaBins*std::log(maxKinEnergy/e)
1094 /std::log(maxKinEnergy/minKinEnergy));
1095 minKinEnergy = e;
1096}
1097
1098//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1099
1101{
1102 nLambdaBins = G4lrint(nLambdaBins*std::log(e/minKinEnergy)
1103 /std::log(maxKinEnergy/minKinEnergy));
1104 maxKinEnergy = e;
1105}
1106
1107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fIsCrossSectionPrim
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
G4ProcessType
@ idxG4ElectronCut
@ idxG4GammaCut
@ idxG4PositronCut
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4bool ForcedInteractionRegion(G4int coupleIdx)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
void UpdateEmModel(const G4String &, G4double, G4double)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4int NumberOfModels() const
void SetFluoFlag(G4bool val)
void DumpModelList(G4int verb)
G4VEmModel * GetModel(G4int, G4bool ver=false)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
const G4DataVector * Initialise(const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4int)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
const std::vector< G4double > * GetDensityFactors()
void InitialiseBaseMaterials(G4PhysicsTable *table)
const std::vector< G4int > * GetCoupleIndexes()
G4bool GetFlag(size_t idx) const
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
G4VAtomDeexcitation * AtomDeexcitation()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
void InitializeForPostStep(const G4Track &)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
size_t length() const
void clearAndDestroy()
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
size_t GetVectorLength() const
G4double GetMaxEnergy() const
G4double Energy(size_t index) const
void FillSecondDerivatives()
void SetSpline(G4bool)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetSafety() const
Definition: G4Step.hh:78
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:613
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:620
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:240
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:186
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:391
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:318
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false)
G4double MeanFreePath(const G4Track &track)
void StartTracking(G4Track *)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
virtual void PrintInfo()=0
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:89
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
void SetMinKinEnergy(G4double e)
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
void SetEmModel(G4VEmModel *, G4int index=1)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4VEmModel * EmModel(G4int index=1)
G4double GetLambda(G4double &kinEnergy, const G4MaterialCutsCouple *couple)
void PrintInfoDefinition()
void UpdateEmModel(const G4String &, G4double, G4double)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void PreparePhysicsTable(const G4ParticleDefinition &)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
G4VEmModel * Model(G4int index=1)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4ParticleChangeForGamma fParticleChange
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
void SetParticle(const G4ParticleDefinition *p)
void SetMaxKinEnergy(G4double e)
const G4Element * GetCurrentElement() const
void SetModel(G4VEmModel *, G4int index=1)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual ~G4VEmProcess()
G4double GetParentWeight() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4TrackStatus GetTrackStatus() const
G4double currentInteractionLength
Definition: G4VProcess.hh:297
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:408
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:214
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
int G4lrint(double ad)
Definition: templates.hh:163
#define DBL_MAX
Definition: templates.hh:83