Geant4_10
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: G4VEmProcess.cc 76957 2013-11-19 15:07:20Z gcosmo $
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"
66 #include "G4PhysicalConstants.hh"
67 #include "G4SystemOfUnits.hh"
68 #include "G4ProcessManager.hh"
69 #include "G4LossTableManager.hh"
70 #include "G4LossTableBuilder.hh"
71 #include "G4Step.hh"
72 #include "G4ParticleDefinition.hh"
73 #include "G4VEmModel.hh"
74 #include "G4DataVector.hh"
75 #include "G4PhysicsTable.hh"
76 #include "G4PhysicsLogVector.hh"
77 #include "G4VParticleChange.hh"
78 #include "G4ProductionCutsTable.hh"
79 #include "G4Region.hh"
80 #include "G4Gamma.hh"
81 #include "G4Electron.hh"
82 #include "G4Positron.hh"
83 #include "G4PhysicsTableHelper.hh"
84 #include "G4EmBiasingManager.hh"
85 #include "G4GenericIon.hh"
86 #include "G4Log.hh"
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
91  G4VDiscreteProcess(name, type),
92  secondaryParticle(0),
93  buildLambdaTable(true),
94  numberOfModels(0),
95  theLambdaTable(0),
96  theLambdaTablePrim(0),
97  theDensityFactor(0),
98  theDensityIdx(0),
99  integral(false),
100  applyCuts(false),
101  startFromNull(false),
102  splineFlag(true),
103  currentModel(0),
104  particle(0),
105  currentParticle(0),
106  currentCouple(0)
107 {
108  SetVerboseLevel(1);
109 
110  // Size of tables assuming spline
111  minKinEnergy = 0.1*keV;
112  maxKinEnergy = 10.0*TeV;
113  nLambdaBins = 77;
114  minKinEnergyPrim = DBL_MAX;
115 
116  // default lambda factor
117  lambdaFactor = 0.8;
118 
119  // default limit on polar angle
120  polarAngleLimit = 0.0;
121  biasFactor = 1.0;
122 
123  // particle types
124  theGamma = G4Gamma::Gamma();
125  theElectron = G4Electron::Electron();
126  thePositron = G4Positron::Positron();
127 
130  secParticles.reserve(5);
131 
132  preStepLambda = 0.0;
133  mfpKinEnergy = DBL_MAX;
134 
135  idxLambda = idxLambdaPrim = 0;
136 
137  modelManager = new G4EmModelManager();
138  biasManager = 0;
139  biasFlag = false;
140  weightFlag = false;
141  lManager = G4LossTableManager::Instance();
142  lManager->Register(this);
143  secID = fluoID = augerID = biasID = -1;
144  mainSecondaries = 100;
145  if("phot" == GetProcessName() || "compt" == GetProcessName()) {
146  mainSecondaries = 1;
147  }
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
153 {
154  if(1 < verboseLevel) {
155  G4cout << "G4VEmProcess destruct " << GetProcessName()
156  << " " << this << " " << theLambdaTable <<G4endl;
157  }
158  if(lManager->IsMaster()) {
159  delete theLambdaTable;
160  delete theLambdaTablePrim;
161  }
162  delete modelManager;
163  delete biasManager;
164  lManager->DeRegister(this);
165  //G4cout << "G4VEmProcess removed " << G4endl;
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
170 void G4VEmProcess::Clear()
171 {
172  currentCouple = 0;
173  preStepLambda = 0.0;
174  mfpKinEnergy = DBL_MAX;
175  idxLambda = idxLambdaPrim = 0;
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
181  const G4Material*)
182 {
183  return 0.0;
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187 
189  const G4Region* region)
190 {
192  modelManager->AddEmModel(order, p, fm, region);
193  if(p) { p->SetParticleChange(pParticleChange); }
194 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197 
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  G4VEmModel* p = 0;
210  if(index >= 0 && index < G4int(emModels.size())) { p = emModels[index]; }
211  return p;
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215 
217  G4double emin, G4double emax)
218 {
219  modelManager->UpdateEmModel(nam, emin, emax);
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
223 
225 {
226  return modelManager->GetModel(idx, ver);
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
230 
232 {
233  G4bool isMaster = false;
234  if(GetMasterProcess() == this) { isMaster = true; }
235  if(!particle) { SetParticle(&part); }
236 
237  if(part.GetParticleType() == "nucleus" &&
238  part.GetParticleSubType() == "generic") {
239 
240  G4String pname = part.GetParticleName();
241  if(pname != "deuteron" && pname != "triton" &&
242  pname != "alpha" && pname != "He3" &&
243  pname != "alpha+" && pname != "helium" &&
244  pname != "hydrogen") {
245 
246  particle = G4GenericIon::GenericIon();
247  }
248  }
249 
250  if(1 < verboseLevel) {
251  G4cout << "G4VEmProcess::PreparePhysicsTable() for "
252  << GetProcessName()
253  << " and particle " << part.GetParticleName()
254  << " local particle " << particle->GetParticleName()
255  << G4endl;
256  }
257 
258  if(particle != &part) { return; }
259 
260  G4LossTableBuilder* bld = lManager->GetTableBuilder();
261 
262  lManager->PreparePhysicsTable(&part, this, isMaster);
263 
264  Clear();
265  InitialiseProcess(particle);
266 
267  const G4ProductionCutsTable* theCoupleTable=
269  size_t n = theCoupleTable->GetTableSize();
270 
271  theEnergyOfCrossSectionMax.resize(n, 0.0);
272  theCrossSectionMax.resize(n, DBL_MAX);
273 
274  // initialisation of models
275  numberOfModels = modelManager->NumberOfModels();
276  for(G4int i=0; i<numberOfModels; ++i) {
277  G4VEmModel* mod = modelManager->GetModel(i);
278  if(0 == i) { currentModel = mod; }
279  mod->SetPolarAngleLimit(polarAngleLimit);
280  mod->SetMasterThread(isMaster);
281  if(mod->HighEnergyLimit() > maxKinEnergy) {
282  mod->SetHighEnergyLimit(maxKinEnergy);
283  }
284  }
285 
286  if(lManager->AtomDeexcitation()) { modelManager->SetFluoFlag(true); }
287  theCuts = modelManager->Initialise(particle,secondaryParticle,
288  2.,verboseLevel);
289  theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
290  theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
291  theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
292 
293  // prepare tables
294  if(buildLambdaTable && isMaster){
295  theLambdaTable =
297  bld->InitialiseBaseMaterials(theLambdaTable);
298  }
299  // high energy table
300  if(isMaster && minKinEnergyPrim < maxKinEnergy){
301  theLambdaTablePrim =
302  G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTablePrim);
303  bld->InitialiseBaseMaterials(theLambdaTablePrim);
304  }
305  // forced biasing
306  if(biasManager) {
307  biasManager->Initialise(part,GetProcessName(),verboseLevel);
308  biasFlag = false;
309  }
310  // defined ID of secondary particles
311  if(isMaster) {
312  G4String nam1 = GetProcessName();
313  G4String nam2 = nam1 + "_fluo" ;
314  G4String nam3 = nam1 + "_auger";
315  G4String nam4 = nam1 + "_split";
316  secID = G4PhysicsModelCatalog::Register(nam1);
317  fluoID = G4PhysicsModelCatalog::Register(nam2);
318  augerID = G4PhysicsModelCatalog::Register(nam3);
319  biasID = G4PhysicsModelCatalog::Register(nam4);
320  }
321 }
322 
323 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324 
326 {
327  const G4VEmProcess* masterProc = 0;
328  if(GetMasterProcess() != this) {
329  masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());
330  }
331 
332  G4String num = part.GetParticleName();
333  if(1 < verboseLevel) {
334  G4cout << "G4VEmProcess::BuildPhysicsTable() for "
335  << GetProcessName()
336  << " and particle " << num
337  << " buildLambdaTable= " << buildLambdaTable
338  << G4endl;
339  }
340 
341  if(particle == &part) {
342 
343  G4LossTableBuilder* bld = lManager->GetTableBuilder();
344 
345  // worker initialisation
346  if(masterProc) {
347  theLambdaTable = masterProc->LambdaTable();
348  theLambdaTablePrim = masterProc->LambdaTablePrim();
349 
350  if(theLambdaTable) {
351  bld->InitialiseBaseMaterials(theLambdaTable);
352  } else if(theLambdaTablePrim) {
353  bld->InitialiseBaseMaterials(theLambdaTablePrim);
354  }
355  theDensityFactor = bld->GetDensityFactors();
356  theDensityIdx = bld->GetCoupleIndexes();
357  if(theLambdaTable) { FindLambdaMax(); }
358 
359  // local initialisation of models
360  G4bool printing = true;
361  numberOfModels = modelManager->NumberOfModels();
362  for(G4int i=0; i<numberOfModels; ++i) {
363  G4VEmModel* mod = GetModelByIndex(i, printing);
364  G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing);
365  mod->InitialiseLocal(particle, mod0);
366  }
367  // master thread
368  } else {
369  theDensityFactor = bld->GetDensityFactors();
370  theDensityIdx = bld->GetCoupleIndexes();
371  if(buildLambdaTable || minKinEnergyPrim < maxKinEnergy) {
372  BuildLambdaTable();
373  }
374  }
375  }
376 
377  // explicitly defined printout by particle name
378  if(1 < verboseLevel ||
379  (0 < verboseLevel && (num == "gamma" || num == "e-" ||
380  num == "e+" || num == "mu+" ||
381  num == "mu-" || num == "proton"||
382  num == "pi+" || num == "pi-" ||
383  num == "kaon+" || num == "kaon-" ||
384  num == "alpha" || num == "anti_proton" ||
385  num == "GenericIon")))
386  {
387  PrintInfoProcess(part);
388  }
389 
390  if(1 < verboseLevel) {
391  G4cout << "G4VEmProcess::BuildPhysicsTable() done for "
392  << GetProcessName()
393  << " and particle " << num
394  << G4endl;
395  }
396 }
397 
398 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
399 
400 void G4VEmProcess::BuildLambdaTable()
401 {
402  if(1 < verboseLevel) {
403  G4cout << "G4EmProcess::BuildLambdaTable() for process "
404  << GetProcessName() << " and particle "
405  << particle->GetParticleName() << " " << this
406  << G4endl;
407  }
408 
409  // Access to materials
410  const G4ProductionCutsTable* theCoupleTable=
412  size_t numOfCouples = theCoupleTable->GetTableSize();
413 
414  G4LossTableBuilder* bld = lManager->GetTableBuilder();
415 
416  G4PhysicsLogVector* aVector = 0;
417  G4PhysicsLogVector* aVectorPrim = 0;
418  G4PhysicsLogVector* bVectorPrim = 0;
419 
420  G4double scale = G4Log(maxKinEnergy/minKinEnergy);
421  G4double emax1 = maxKinEnergy;
422  if(minKinEnergyPrim < maxKinEnergy) { emax1 = minKinEnergyPrim; }
423 
424  for(size_t i=0; i<numOfCouples; ++i) {
425 
426  if (bld->GetFlag(i)) {
427 
428  // create physics vector and fill it
429  const G4MaterialCutsCouple* couple =
430  theCoupleTable->GetMaterialCutsCouple(i);
431 
432  // build main table
433  if(buildLambdaTable) {
434  delete (*theLambdaTable)[i];
435 
436  // if start from zero then change the scale
437  G4double emin = minKinEnergy;
438  G4bool startNull = false;
439  if(startFromNull) {
440  G4double e = MinPrimaryEnergy(particle,couple->GetMaterial());
441  if(e >= emin) {
442  emin = e;
443  startNull = true;
444  }
445  }
446  G4double emax = emax1;
447  if(emax <= emin) { emax = 2*emin; }
448  G4int bin = G4lrint(nLambdaBins*G4Log(emax/emin)/scale);
449  if(bin < 3) { bin = 3; }
450  aVector = new G4PhysicsLogVector(emin, emax, bin);
451  aVector->SetSpline(splineFlag);
452  modelManager->FillLambdaVector(aVector, couple, startNull);
453  if(splineFlag) { aVector->FillSecondDerivatives(); }
454  G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector);
455  }
456  // build high energy table
457  if(minKinEnergyPrim < maxKinEnergy) {
458  delete (*theLambdaTablePrim)[i];
459 
460  // start not from zero
461  if(!bVectorPrim) {
462  G4int bin =
463  G4lrint(nLambdaBins*G4Log(maxKinEnergy/minKinEnergyPrim)/scale);
464  if(bin < 3) { bin = 3; }
465  aVectorPrim =
466  new G4PhysicsLogVector(minKinEnergyPrim, maxKinEnergy, bin);
467  bVectorPrim = aVectorPrim;
468  } else {
469  aVectorPrim = new G4PhysicsLogVector(*bVectorPrim);
470  }
471  // always use spline
472  aVectorPrim->SetSpline(true);
473  modelManager->FillLambdaVector(aVectorPrim, couple, false,
475  aVectorPrim->FillSecondDerivatives();
476  G4PhysicsTableHelper::SetPhysicsVector(theLambdaTablePrim, i,
477  aVectorPrim);
478  }
479  }
480  }
481 
482  if(buildLambdaTable) { FindLambdaMax(); }
483 
484  if(1 < verboseLevel) {
485  G4cout << "Lambda table is built for "
486  << particle->GetParticleName()
487  << G4endl;
488  }
489 }
490 
491 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492 
493 void G4VEmProcess::PrintInfoProcess(const G4ParticleDefinition& part)
494 {
495  if(verboseLevel > 0) {
496  G4cout << G4endl << GetProcessName() << ": for "
497  << part.GetParticleName();
498  if(integral) { G4cout << ", integral: 1 "; }
499  if(applyCuts) { G4cout << ", applyCuts: 1 "; }
500  G4cout << " SubType= " << GetProcessSubType();;
501  if(biasFactor != 1.0) { G4cout << " BiasingFactor= " << biasFactor; }
502  G4cout << G4endl;
503  if(buildLambdaTable) {
504  size_t length = theLambdaTable->length();
505  for(size_t i=0; i<length; ++i) {
506  G4PhysicsVector* v = (*theLambdaTable)[i];
507  if(v) {
508  G4cout << " Lambda table from "
509  << G4BestUnit(minKinEnergy,"Energy")
510  << " to "
511  << G4BestUnit(v->GetMaxEnergy(),"Energy")
512  << " in " << v->GetVectorLength()-1
513  << " bins, spline: "
514  << splineFlag
515  << G4endl;
516  break;
517  }
518  }
519  }
520  if(minKinEnergyPrim < maxKinEnergy) {
521  size_t length = theLambdaTablePrim->length();
522  for(size_t i=0; i<length; ++i) {
523  G4PhysicsVector* v = (*theLambdaTablePrim)[i];
524  if(v) {
525  G4cout << " LambdaPrime table from "
526  << G4BestUnit(minKinEnergyPrim,"Energy")
527  << " to "
528  << G4BestUnit(maxKinEnergy,"Energy")
529  << " in " << v->GetVectorLength()-1
530  << " bins "
531  << G4endl;
532  break;
533  }
534  }
535  }
536  PrintInfo();
537  modelManager->DumpModelList(verboseLevel);
538  }
539 
540  if(verboseLevel > 2 && buildLambdaTable) {
541  G4cout << " LambdaTable address= " << theLambdaTable << G4endl;
542  if(theLambdaTable) { G4cout << (*theLambdaTable) << G4endl; }
543  }
544 }
545 
546 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547 
549 {
550  // reset parameters for the new track
551  currentParticle = track->GetParticleDefinition();
553  //currentInteractionLength = -1.0;
554  // theInitialNumberOfInteractionLength=-1.0;
555  mfpKinEnergy = DBL_MAX;
556 
557  // forced biasing only for primary particles
558  if(biasManager) {
559  if(0 == track->GetParentID()) {
560  // primary particle
561  biasFlag = true;
562  biasManager->ResetForcedInteraction();
563  }
564  }
565 }
566 
567 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
568 
570  const G4Track& track,
571  G4double previousStepSize,
573 {
574  *condition = NotForced;
575  G4double x = DBL_MAX;
576 
577  preStepKinEnergy = track.GetKineticEnergy();
578  DefineMaterial(track.GetMaterialCutsCouple());
579  SelectModel(preStepKinEnergy, currentCoupleIndex);
580 
581  if(!currentModel->IsActive(preStepKinEnergy)) {
583  return x;
584  }
585 
586  // forced biasing only for primary particles
587  if(biasManager) {
588  if(0 == track.GetParentID()) {
589  if(biasFlag &&
590  biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
591  return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
592  }
593  }
594  }
595 
596  // compute mean free path
597  if(preStepKinEnergy < mfpKinEnergy) {
598  if (integral) { ComputeIntegralLambda(preStepKinEnergy); }
599  else { preStepLambda = GetCurrentLambda(preStepKinEnergy); }
600 
601  // zero cross section
602  if(preStepLambda <= 0.0) {
605  }
606  }
607 
608  // non-zero cross section
609  if(preStepLambda > 0.0) {
610 
612 
613  // beggining of tracking (or just after DoIt of this process)
614  // ResetNumberOfInteractionLengthLeft();
615 
618 
619  } else if(currentInteractionLength < DBL_MAX) {
620 
621  // subtract NumberOfInteractionLengthLeft using previous step
623  previousStepSize/currentInteractionLength;
624  //SubtractNumberOfInteractionLengthLeft(previousStepSize);
627  }
628  }
629 
630  // new mean free path and step limit for the next step
631  currentInteractionLength = 1.0/preStepLambda;
633 #ifdef G4VERBOSE
634  if (verboseLevel>2){
635  G4cout << "G4VEmProcess::PostStepGetPhysicalInteractionLength ";
636  G4cout << "[ " << GetProcessName() << "]" << G4endl;
637  G4cout << " for " << currentParticle->GetParticleName()
638  << " in Material " << currentMaterial->GetName()
639  << " Ekin(MeV)= " << preStepKinEnergy/MeV
640  <<G4endl;
641  G4cout << " MeanFreePath = " << currentInteractionLength/cm << "[cm]"
642  << " InteractionLength= " << x/cm <<"[cm] " <<G4endl;
643  }
644 #endif
645  }
646  return x;
647 }
648 
649 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
650 
652  const G4Step& step)
653 {
654  // In all cases clear number of interaction lengths
656  mfpKinEnergy = DBL_MAX;
657 
659 
660  // Do not make anything if particle is stopped, the annihilation then
661  // should be performed by the AtRestDoIt!
662  if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
663 
664  G4double finalT = track.GetKineticEnergy();
665 
666  // forced process - should happen only once per track
667  if(biasFlag) {
668  if(biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
669  biasFlag = false;
670  }
671  }
672 
673  // Integral approach
674  if (integral) {
675  G4double lx = GetLambda(finalT, currentCouple);
676  if(preStepLambda<lx && 1 < verboseLevel) {
677  G4cout << "WARNING: for " << currentParticle->GetParticleName()
678  << " and " << GetProcessName()
679  << " E(MeV)= " << finalT/MeV
680  << " preLambda= " << preStepLambda << " < "
681  << lx << " (postLambda) "
682  << G4endl;
683  }
684 
685  if(preStepLambda*G4UniformRand() > lx) {
687  return &fParticleChange;
688  }
689  }
690 
691  SelectModel(finalT, currentCoupleIndex);
692  if(!currentModel->IsActive(finalT)) { return &fParticleChange; }
693 
694  // define new weight for primary and secondaries
696  if(weightFlag) {
697  weight /= biasFactor;
699  }
700 
701  /*
702  if(0 < verboseLevel) {
703  G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
704  << finalT/MeV
705  << " MeV; model= (" << currentModel->LowEnergyLimit()
706  << ", " << currentModel->HighEnergyLimit() << ")"
707  << G4endl;
708  }
709  */
710 
711  // sample secondaries
712  secParticles.clear();
713  currentModel->SampleSecondaries(&secParticles,
714  currentCouple,
715  track.GetDynamicParticle(),
716  (*theCuts)[currentCoupleIndex]);
717 
718  G4int num0 = secParticles.size();
719 
720  // splitting or Russian roulette
721  if(biasManager) {
722  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
723  G4double eloss = 0.0;
724  weight *= biasManager->ApplySecondaryBiasing(
725  secParticles, track, currentModel, &fParticleChange, eloss,
726  currentCoupleIndex, (*theCuts)[currentCoupleIndex],
727  step.GetPostStepPoint()->GetSafety());
728  if(eloss > 0.0) {
731  }
732  }
733  }
734 
735  // save secondaries
736  G4int num = secParticles.size();
737  if(num > 0) {
738 
741  G4double time = track.GetGlobalTime();
742 
743  for (G4int i=0; i<num; ++i) {
744  if (secParticles[i]) {
745  G4DynamicParticle* dp = secParticles[i];
747  G4double e = dp->GetKineticEnergy();
748  G4bool good = true;
749  if(applyCuts) {
750  if (p == theGamma) {
751  if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
752 
753  } else if (p == theElectron) {
754  if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
755 
756  } else if (p == thePositron) {
757  if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
758  e < (*theCutsPositron)[currentCoupleIndex]) {
759  good = false;
760  e += 2.0*electron_mass_c2;
761  }
762  }
763  // added secondary if it is good
764  }
765  if (good) {
766  G4Track* t = new G4Track(dp, time, track.GetPosition());
768  t->SetWeight(weight);
770 
771  // define type of secondary
772  if(i < mainSecondaries) { t->SetCreatorModelIndex(secID); }
773  else if(i < num0) {
774  if(p == theGamma) {
775  t->SetCreatorModelIndex(fluoID);
776  } else {
777  t->SetCreatorModelIndex(augerID);
778  }
779  } else {
780  t->SetCreatorModelIndex(biasID);
781  }
782 
783  //G4cout << "Secondary(post step) has weight " << t->GetWeight()
784  // << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
785  } else {
786  delete dp;
787  edep += e;
788  }
789  }
790  }
792  }
793 
796  if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
799  }
800 
801  // ClearNumberOfInteractionLengthLeft();
802  return &fParticleChange;
803 }
804 
805 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
806 
808  const G4String& directory,
809  G4bool ascii)
810 {
811  G4bool yes = true;
812 
813  if ( theLambdaTable && part == particle) {
814  const G4String name =
815  GetPhysicsTableFileName(part,directory,"Lambda",ascii);
816  yes = theLambdaTable->StorePhysicsTable(name,ascii);
817 
818  if ( yes ) {
819  G4cout << "Physics table is stored for " << particle->GetParticleName()
820  << " and process " << GetProcessName()
821  << " in the directory <" << directory
822  << "> " << G4endl;
823  } else {
824  G4cout << "Fail to store Physics Table for "
825  << particle->GetParticleName()
826  << " and process " << GetProcessName()
827  << " in the directory <" << directory
828  << "> " << G4endl;
829  }
830  }
831  if ( theLambdaTablePrim && part == particle) {
832  const G4String name =
833  GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
834  yes = theLambdaTablePrim->StorePhysicsTable(name,ascii);
835 
836  if ( yes ) {
837  G4cout << "Physics table prim is stored for "
838  << particle->GetParticleName()
839  << " and process " << GetProcessName()
840  << " in the directory <" << directory
841  << "> " << G4endl;
842  } else {
843  G4cout << "Fail to store Physics Table Prim for "
844  << particle->GetParticleName()
845  << " and process " << GetProcessName()
846  << " in the directory <" << directory
847  << "> " << G4endl;
848  }
849  }
850  return yes;
851 }
852 
853 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
854 
856  const G4String& directory,
857  G4bool ascii)
858 {
859  if(1 < verboseLevel) {
860  G4cout << "G4VEmProcess::RetrievePhysicsTable() for "
861  << part->GetParticleName() << " and process "
862  << GetProcessName() << G4endl;
863  }
864  G4bool yes = true;
865 
866  if((!buildLambdaTable && minKinEnergyPrim > maxKinEnergy)
867  || particle != part) { return yes; }
868 
869  const G4String particleName = part->GetParticleName();
870  G4String filename;
871 
872  if(buildLambdaTable) {
873  filename = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
874  yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTable,
875  filename,ascii);
876  if ( yes ) {
877  if (0 < verboseLevel) {
878  G4cout << "Lambda table for " << particleName
879  << " is Retrieved from <"
880  << filename << ">"
881  << G4endl;
882  }
883  if(lManager->SplineFlag()) {
884  size_t n = theLambdaTable->length();
885  for(size_t i=0; i<n; ++i) {
886  if((* theLambdaTable)[i]) {
887  (* theLambdaTable)[i]->SetSpline(true);
888  }
889  }
890  }
891  } else {
892  if (1 < verboseLevel) {
893  G4cout << "Lambda table for " << particleName << " in file <"
894  << filename << "> is not exist"
895  << G4endl;
896  }
897  }
898  }
899  if(minKinEnergyPrim < maxKinEnergy) {
900  filename = GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
901  yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTablePrim,
902  filename,ascii);
903  if ( yes ) {
904  if (0 < verboseLevel) {
905  G4cout << "Lambda table prim for " << particleName
906  << " is Retrieved from <"
907  << filename << ">"
908  << G4endl;
909  }
910  if(lManager->SplineFlag()) {
911  size_t n = theLambdaTablePrim->length();
912  for(size_t i=0; i<n; ++i) {
913  if((* theLambdaTablePrim)[i]) {
914  (* theLambdaTablePrim)[i]->SetSpline(true);
915  }
916  }
917  }
918  } else {
919  if (1 < verboseLevel) {
920  G4cout << "Lambda table prim for " << particleName << " in file <"
921  << filename << "> is not exist"
922  << G4endl;
923  }
924  }
925  }
926 
927  return yes;
928 }
929 
930 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
931 
932 G4double
934  const G4MaterialCutsCouple* couple)
935 {
936  // Cross section per atom is calculated
937  DefineMaterial(couple);
938  G4double cross = 0.0;
939  if(buildLambdaTable && theLambdaTable) {
940  cross = GetCurrentLambda(kineticEnergy);
941  } else {
942  SelectModel(kineticEnergy, currentCoupleIndex);
943  cross = fFactor*currentModel->CrossSectionPerVolume(currentMaterial,
944  currentParticle,
945  kineticEnergy);
946  }
947 
948  if(cross < 0.0) { cross = 0.0; }
949  return cross;
950 }
951 
952 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
953 
955  G4double,
957 {
958  *condition = NotForced;
959  return G4VEmProcess::MeanFreePath(track);
960 }
961 
962 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
963 
965 {
966  DefineMaterial(track.GetMaterialCutsCouple());
967  preStepLambda = GetCurrentLambda(track.GetKineticEnergy());
968  G4double x = DBL_MAX;
969  if(0.0 < preStepLambda) { x = 1.0/preStepLambda; }
970  return x;
971 }
972 
973 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
974 
975 G4double
977  G4double Z, G4double A, G4double cut)
978 {
979  SelectModel(kineticEnergy, currentCoupleIndex);
980  G4double x = 0.0;
981  if(currentModel) {
982  x = currentModel->ComputeCrossSectionPerAtom(currentParticle,kineticEnergy,
983  Z,A,cut);
984  }
985  return x;
986 }
987 
988 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
989 
990 void G4VEmProcess::FindLambdaMax()
991 {
992  if(1 < verboseLevel) {
993  G4cout << "### G4VEmProcess::FindLambdaMax: "
994  << particle->GetParticleName()
995  << " and process " << GetProcessName() << " " << G4endl;
996  }
997  size_t n = theLambdaTable->length();
998  G4PhysicsVector* pv;
999  G4double e, ss, emax, smax;
1000 
1001  size_t i;
1002 
1003  // first loop on existing vectors
1004  for (i=0; i<n; ++i) {
1005  pv = (*theLambdaTable)[i];
1006  if(pv) {
1007  size_t nb = pv->GetVectorLength();
1008  emax = DBL_MAX;
1009  smax = 0.0;
1010  if(nb > 0) {
1011  for (size_t j=0; j<nb; ++j) {
1012  e = pv->Energy(j);
1013  ss = (*pv)(j);
1014  if(ss > smax) {
1015  smax = ss;
1016  emax = e;
1017  }
1018  }
1019  }
1020  theEnergyOfCrossSectionMax[i] = emax;
1021  theCrossSectionMax[i] = smax;
1022  if(1 < verboseLevel) {
1023  G4cout << "For " << particle->GetParticleName()
1024  << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
1025  << " lambda= " << smax << G4endl;
1026  }
1027  }
1028  }
1029  // second loop using base materials
1030  for (i=0; i<n; ++i) {
1031  pv = (*theLambdaTable)[i];
1032  if(!pv){
1033  G4int j = (*theDensityIdx)[i];
1034  theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
1035  theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
1036  }
1037  }
1038 }
1039 
1040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1041 
1044 {
1045  G4PhysicsVector* v =
1046  new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nLambdaBins);
1047  v->SetSpline(lManager->SplineFlag());
1048  return v;
1049 }
1050 
1051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1052 
1054 {
1055  const G4Element* elm = 0;
1056  if(currentModel) {elm = currentModel->GetCurrentElement(); }
1057  return elm;
1058 }
1059 
1060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1061 
1063 {
1064  if(f > 0.0) {
1065  biasFactor = f;
1066  weightFlag = flag;
1067  if(1 < verboseLevel) {
1068  G4cout << "### SetCrossSectionBiasingFactor: for "
1069  << particle->GetParticleName()
1070  << " and process " << GetProcessName()
1071  << " biasFactor= " << f << " weightFlag= " << flag
1072  << G4endl;
1073  }
1074  }
1075 }
1076 
1077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1078 
1079 void
1081  G4bool flag)
1082 {
1083  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1084  if(1 < verboseLevel) {
1085  G4cout << "### ActivateForcedInteraction: for "
1086  << particle->GetParticleName()
1087  << " and process " << GetProcessName()
1088  << " length(mm)= " << length/mm
1089  << " in G4Region <" << r
1090  << "> weightFlag= " << flag
1091  << G4endl;
1092  }
1093  weightFlag = flag;
1094  biasManager->ActivateForcedInteraction(length, r);
1095 }
1096 
1097 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1098 
1099 void
1101  G4double factor,
1102  G4double energyLimit)
1103 {
1104  if (0.0 <= factor) {
1105 
1106  // Range cut can be applied only for e-
1107  if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1108  { return; }
1109 
1110  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1111  biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1112  if(1 < verboseLevel) {
1113  G4cout << "### ActivateSecondaryBiasing: for "
1114  << " process " << GetProcessName()
1115  << " factor= " << factor
1116  << " in G4Region <" << region
1117  << "> energyLimit(MeV)= " << energyLimit/MeV
1118  << G4endl;
1119  }
1120  }
1121 }
1122 
1123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1124 
1126 {
1127  nLambdaBins = G4lrint(nLambdaBins*std::log(maxKinEnergy/e)
1128  /std::log(maxKinEnergy/minKinEnergy));
1129  minKinEnergy = e;
1130 }
1131 
1132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1133 
1135 {
1136  nLambdaBins = G4lrint(nLambdaBins*std::log(e/minKinEnergy)
1137  /std::log(maxKinEnergy/minKinEnergy));
1138  maxKinEnergy = e;
1139 }
1140 
1141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
G4double condition(const G4ErrorSymMatrix &m)
G4PhysicsTable * LambdaTable() const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
virtual void PrintInfo()=0
G4int GetParentID() const
G4bool SplineFlag() const
G4double GetLambda(G4double &kinEnergy, const G4MaterialCutsCouple *couple)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4PhysicsTable * LambdaTablePrim() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
const std::vector< G4double > * GetDensityFactors()
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
static G4LossTableManager * Instance()
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4double GetMaxEnergy() const
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4int verboseLevel
Definition: G4VProcess.hh:368
void UpdateEmModel(const G4String &, G4double, G4double)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetKineticEnergy() const
Int_t index
Definition: macro.C:9
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:90
const G4DynamicParticle * GetDynamicParticle() const
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4VEmModel * EmModel(G4int index=1) const
G4VEmModel * GetModel(G4int, G4bool ver=false)
void DeRegister(G4VEnergyLossProcess *p)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:176
const G4ThreeVector & GetPosition() const
G4ParticleChangeForGamma fParticleChange
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:206
G4TrackStatus GetTrackStatus() const
G4bool ForcedInteractionRegion(G4int coupleIdx)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
float bin[41]
Definition: plottest35.C:14
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
virtual ~G4VEmProcess()
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4bool GetFlag(size_t idx) const
const XML_Char * name
Definition: expat.h:151
void SetTouchableHandle(const G4TouchableHandle &apValue)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4double GetParentWeight() const
tuple x
Definition: test.py:50
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
size_t GetVectorLength() const
const G4String & GetParticleSubType() const
double weight
Definition: plottest35.C:25
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
Double_t edep
Definition: macro.C:13
void UpdateEmModel(const G4String &, G4double, G4double)
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
const G4int smax
const G4String & GetParticleName() const
TFile f
Definition: plotHisto.C:6
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4LossTableBuilder * GetTableBuilder()
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetCreatorModelIndex(G4int idx)
Definition: G4Track.hh:245
void SetSpline(G4bool)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
G4double MeanFreePath(const G4Track &track)
void ProposeWeight(G4double finalWeight)
void SetSecondaryWeightByProcess(G4bool)
Char_t n[5]
void SetEmModel(G4VEmModel *, G4int index=1)
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
Float_t Z
Definition: plot.C:39
bool G4bool
Definition: G4Types.hh:79
void SetParticle(const G4ParticleDefinition *p)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
TString part[npart]
Definition: Style.C:32
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
const G4ParticleDefinition * GetParticleDefinition() const
G4bool IsMaster() const
const G4String & GetParticleType() const
void Register(G4VEnergyLossProcess *p)
float electron_mass_c2
Definition: hepunit.py:274
Definition: G4Step.hh:76
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetGlobalTime() const
void PreparePhysicsTable(const G4ParticleDefinition &)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double Energy(size_t index) const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
Double_t scale
Definition: plot.C:11
void AddSecondary(G4Track *aSecondary)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
tuple v
Definition: test.py:18
const G4TouchableHandle & GetTouchableHandle() const
string pname
Definition: eplot.py:33
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
void DumpModelList(G4int verb)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:669
size_t length() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4int size() const
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
static G4ProductionCutsTable * GetProductionCutsTable()
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
#define fm
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
jump r
Definition: plot.C:36
G4int NumberOfModels() const
void SetMaxKinEnergy(G4double e)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetProposedKineticEnergy() const
void InitialiseBaseMaterials(G4PhysicsTable *table)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double GetLocalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4Element * GetCurrentElement() const
G4double GetSafety() const
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void SetFluoFlag(G4bool val)
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4TrackStatus GetTrackStatus() const
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4VAtomDeexcitation * AtomDeexcitation()
void SetMinKinEnergy(G4double e)
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:300
static G4int Register(G4String &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void InitializeForPostStep(const G4Track &)
G4ForceCondition
void StartTracking(G4Track *)
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
#define DBL_MAX
Definition: templates.hh:83
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
const G4Material * GetMaterial() const
const std::vector< G4int > * GetCoupleIndexes()
G4ProcessType
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440