Geant4  10.02.p03
G4AdjointPhysicsList.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 //
28 //
29 // $Id: G4AdjointPhysicsList.cc 102356 2017-01-23 16:22:42Z gcosmo $
30 //
32 // Class Name: G4AdjointPhysicsList
33 // Author: L. Desorgher
34 // Organisation: SpaceIT GmbH
35 // Contract: ESA contract 21435/08/NL/AT
36 // Customer: ESA/ESTEC
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41 
42 #include "G4AdjointPhysicsList.hh"
43 #include "G4ProcessManager.hh"
44 #include "G4ParticleTypes.hh"
46 #include "G4SystemOfUnits.hh"
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
52  fEminusIonisation(0),fPIonisation(0),
53  fUse_eionisation(true),fUse_pionisation(true),
54  fUse_brem(true),fUse_compton(true),fUse_ms(true),
55  fUse_egain_fluctuation(true),fUse_peeffect(true),
56  fEmin_adj_models(1.*keV), fEmax_adj_models(1.*MeV),
57  fCS_biasing_factor_compton(1.),fCS_biasing_factor_brem(1.),
58  fCS_biasing_factor_ionisation(1.),fCS_biasing_factor_PEeffect(1.)
59 {
60  defaultCutValue = 1.0*mm;
61  SetVerboseLevel(1);
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
68 {
69 }
71 {
72  // In this method, static member functions should be called
73  // for all particles which you want to use.
74  // This ensures that objects of these particle types will be
75  // created in the program.
81 }
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84 
86 {
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {
94  // pseudo-particles
97 
98  // gamma
100 
101  // optical photon
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108 {
109  // leptons
114 
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124 {
125 // mesons
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142 {
143 // barions
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 
152 #include"G4AdjointGamma.hh"
153 #include"G4AdjointElectron.hh"
154 #include"G4AdjointProton.hh"
156 {
157 // adjoint_gammma
159 
160 // adjoint_electron
162 
163 // adjoint_proton
165 }
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
168 
170 {
172  ConstructEM();
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
178 //#include "G4PEEffectFluoModel.hh"
179 #include "G4ComptonScattering.hh"
180 #include "G4GammaConversion.hh"
181 #include "G4PhotoElectricEffect.hh"
182 #include "G4eMultipleScattering.hh"
183 #include "G4hMultipleScattering.hh"
184 #include "G4eIonisation.hh"
185 #include "G4eBremsstrahlung.hh"
186 #include "G4eplusAnnihilation.hh"
187 #include "G4hIonisation.hh"
188 #include "G4ionIonisation.hh"
189 //#include "G4IonParametrisedLossModel.hh"
190 
192 #include "G4eInverseIonisation.hh"
194 #include "G4AdjointCSManager.hh"
197 #include "G4AdjointComptonModel.hh"
198 #include "G4eInverseCompton.hh"
199 #include "G4InversePEEffect.hh"
202 #include "G4hInverseIonisation.hh"
205 #include "G4IonInverseIonisation.hh"
207 
208 #include "G4AdjointSimManager.hh"
210 
212 { G4AdjointCSManager* theCSManager =
214 
215  G4AdjointSimManager* theAdjointSimManager =
217 
218  theCSManager->RegisterAdjointParticle(
220 
222  theCSManager->RegisterAdjointParticle(
224 
225  if (fUse_eionisation) {
227  new G4eIonisation();
229  }
230 
231  if (fUse_pionisation) {
234  theCSManager->RegisterAdjointParticle(
236  }
237 
238  G4eBremsstrahlung* theeminusBremsstrahlung = 0;
240  theeminusBremsstrahlung = new G4eBremsstrahlung();
241 
242  G4ComptonScattering* theComptonScattering =0;
243  if (fUse_compton) theComptonScattering = new G4ComptonScattering();
244 
245  G4PhotoElectricEffect* thePEEffect =0;
246  if (fUse_peeffect) thePEEffect = new G4PhotoElectricEffect();
247 
248  G4eMultipleScattering* theeminusMS = 0;
249  G4hMultipleScattering* thepMS= 0;
250  if (fUse_ms) {
251  theeminusMS = new G4eMultipleScattering();
252  thepMS = new G4hMultipleScattering();
253  }
254 
255  G4VProcess* theGammaConversion =0;
256  if (fUse_gamma_conversion) theGammaConversion = new G4GammaConversion();
257 
258  //Define adjoint e- ionisation
259  //-------------------
260  G4AdjointeIonisationModel* theeInverseIonisationModel = 0;
261  G4eInverseIonisation* theeInverseIonisationProjToProjCase = 0 ;
262  G4eInverseIonisation* theeInverseIonisationProdToProjCase = 0;
263  if (fUse_eionisation) {
264  theeInverseIonisationModel = new G4AdjointeIonisationModel();
265  theeInverseIonisationModel->SetHighEnergyLimit(
267  theeInverseIonisationModel->SetLowEnergyLimit(
269  theeInverseIonisationModel->SetCSBiasingFactor(
271  theeInverseIonisationProjToProjCase =
272  new G4eInverseIonisation(true,"Inv_eIon",
273  theeInverseIonisationModel);
274  theeInverseIonisationProdToProjCase =
275  new G4eInverseIonisation(false,"Inv_eIon1",
276  theeInverseIonisationModel);
277  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
278  }
279 
280  //Define adjoint Bremsstrahlung
281  //-------------------------------
282  G4AdjointBremsstrahlungModel* theeInverseBremsstrahlungModel = 0;
283  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProjToProjCase = 0;
284  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProdToProjCase = 0;
285 
286  if (fUse_brem && fUse_eionisation) {
287  theeInverseBremsstrahlungModel = new G4AdjointBremsstrahlungModel();
288  theeInverseBremsstrahlungModel->SetHighEnergyLimit(fEmax_adj_models);
289  theeInverseBremsstrahlungModel->SetLowEnergyLimit(fEmin_adj_models);
290  theeInverseBremsstrahlungModel->SetCSBiasingFactor(
292  theeInverseBremsstrahlungProjToProjCase = new G4eInverseBremsstrahlung(
293  true,"Inv_eBrem",theeInverseBremsstrahlungModel);
294  theeInverseBremsstrahlungProdToProjCase = new G4eInverseBremsstrahlung(
295  false,"Inv_eBrem1",theeInverseBremsstrahlungModel);
296  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
297  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
298  }
299 
300 
301  //Define adjoint Compton
302  //---------------------
303 
304  G4AdjointComptonModel* theeInverseComptonModel = 0;
305  G4eInverseCompton* theeInverseComptonProjToProjCase = 0;
306  G4eInverseCompton* theeInverseComptonProdToProjCase = 0;
307 
308  if (fUse_compton) {
309  theeInverseComptonModel = new G4AdjointComptonModel();
310  theeInverseComptonModel->SetHighEnergyLimit(fEmax_adj_models);
311  theeInverseComptonModel->SetLowEnergyLimit(fEmin_adj_models);
312  theeInverseComptonModel->SetDirectProcess(theComptonScattering);
313  theeInverseComptonModel->SetUseMatrix(false);
314  theeInverseComptonModel->SetCSBiasingFactor( fCS_biasing_factor_compton);
315  theeInverseComptonProjToProjCase = new G4eInverseCompton(true,"Inv_Compt",
316  theeInverseComptonModel);
317  theeInverseComptonProdToProjCase = new G4eInverseCompton(false,"Inv_Compt1",
318  theeInverseComptonModel);
319  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
320  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
321  }
322 
323  //Define adjoint PEEffect
324  //---------------------
325  G4AdjointPhotoElectricModel* theInversePhotoElectricModel = 0;
326  G4InversePEEffect* theInversePhotoElectricProcess = 0;
327 
328  if (fUse_peeffect) {
329  theInversePhotoElectricModel = new G4AdjointPhotoElectricModel();
330  theInversePhotoElectricModel->SetHighEnergyLimit(fEmax_adj_models);
331  theInversePhotoElectricModel->SetLowEnergyLimit(fEmin_adj_models);
332  theInversePhotoElectricModel->SetCSBiasingFactor(
334  theInversePhotoElectricProcess = new G4InversePEEffect("Inv_PEEffect",
335  theInversePhotoElectricModel);
336  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
337  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
338  }
339 
340 
341  //Define adjoint ionisation for protons
342  //---------------------
343  G4AdjointhIonisationModel* thepInverseIonisationModel = 0;
344  G4hInverseIonisation* thepInverseIonisationProjToProjCase = 0 ;
345  G4hInverseIonisation* thepInverseIonisationProdToProjCase = 0;
346  if (fUse_pionisation) {
347  thepInverseIonisationModel = new G4AdjointhIonisationModel(
348  G4Proton::Proton());
349  thepInverseIonisationModel->SetHighEnergyLimit(fEmax_adj_models);
350  thepInverseIonisationModel->SetLowEnergyLimit(fEmin_adj_models);
351  thepInverseIonisationModel->SetUseMatrix(false);
352  thepInverseIonisationProjToProjCase = new G4hInverseIonisation(true,
353  "Inv_pIon",thepInverseIonisationModel);
354  thepInverseIonisationProdToProjCase = new G4hInverseIonisation(false,
355  "Inv_pIon1",thepInverseIonisationModel);
356  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
357  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("proton"));
358  }
359 
360 
361 
362  //Declare the processes active for the different particles
363  //--------------------------------------------------------
365  particleIterator->reset();
366  while( (*particleIterator)() ){
367  G4ParticleDefinition* particle = particleIterator->value();
368  G4ProcessManager* pmanager = particle->GetProcessManager();
369  if (!pmanager) {
370  pmanager = new G4ProcessManager(particle);
371  particle->SetProcessManager(pmanager);
372  }
373 
374  G4String particleName = particle->GetParticleName();
375  if (particleName == "e-") {
376  if (fUse_ms && fUse_eionisation) pmanager->AddProcess(theeminusMS);
377  if (fUse_eionisation){
378  pmanager->AddProcess(fEminusIonisation);
380  RegisterEnergyLossProcess(fEminusIonisation,particle);
381  }
382  if (fUse_brem && fUse_eionisation) {
383  pmanager->AddProcess(theeminusBremsstrahlung);
385  RegisterEnergyLossProcess(theeminusBremsstrahlung,particle);
386  }
387  G4int n_order=0;
388  if (fUse_ms && fUse_eionisation) {
389  n_order++;
390  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
391  }
392  if (fUse_eionisation) {
393  n_order++;
395  }
396  if (fUse_brem && fUse_eionisation) {
397  n_order++;
398  pmanager->SetProcessOrdering(theeminusBremsstrahlung,
399  idxAlongStep,n_order);
400  }
401  n_order=0;
402  if (fUse_ms && fUse_eionisation) {
403  n_order++;
404  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
405  }
406  if (fUse_eionisation) {
407  n_order++;
409  }
410  if (fUse_brem && fUse_eionisation) {
411  n_order++;
412  pmanager->SetProcessOrdering(theeminusBremsstrahlung,idxPostStep,
413  n_order);
414  }
415  }
416 
417  if (particleName == "adj_e-") {
418  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
419  if (fUse_eionisation ) {
420  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
421  theContinuousGainOfEnergy->SetLossFluctuations(
423  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(
425  theContinuousGainOfEnergy->SetDirectParticle(G4Electron::Electron());
426  pmanager->AddProcess(theContinuousGainOfEnergy);
427  }
428  G4int n_order=0;
429  if (fUse_ms) {
430  n_order++;
431  pmanager->AddProcess(theeminusMS);
432  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
433  }
434 
435  n_order++;
436  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
437  n_order);
438 
439 
440  n_order++;
441  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
443  pmanager->AddProcess(theAlongStepWeightCorrection);
444  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
445  idxAlongStep,
446  n_order);
447  n_order=0;
448  if (fUse_eionisation) {
449  pmanager->AddProcess(theeInverseIonisationProjToProjCase);
450  pmanager->AddProcess(theeInverseIonisationProdToProjCase);
451  n_order++;
452  pmanager->SetProcessOrdering(theeInverseIonisationProjToProjCase,
453  idxPostStep,n_order);
454  n_order++;
455  pmanager->SetProcessOrdering(theeInverseIonisationProdToProjCase,
456  idxPostStep,n_order);
457  }
458 
459  if (fUse_brem && fUse_eionisation) {
460  pmanager->AddProcess(theeInverseBremsstrahlungProjToProjCase);
461  n_order++;
462  pmanager->SetProcessOrdering(
463  theeInverseBremsstrahlungProjToProjCase,
464  idxPostStep,n_order);
465  }
466 
467  if (fUse_compton) {
468  pmanager->AddProcess(theeInverseComptonProdToProjCase);
469  n_order++;
470  pmanager->SetProcessOrdering(theeInverseComptonProdToProjCase,
471  idxPostStep,n_order);
472  }
473  if (fUse_peeffect) {
474  pmanager->AddDiscreteProcess(theInversePhotoElectricProcess);
475  n_order++;
476  pmanager->SetProcessOrdering(theInversePhotoElectricProcess,
477  idxPostStep,n_order);
478  }
479  if (fUse_pionisation) {
480  pmanager->AddProcess(thepInverseIonisationProdToProjCase);
481  n_order++;
482  pmanager->SetProcessOrdering(thepInverseIonisationProdToProjCase,
483  idxPostStep,n_order);
484  }
485  if (fUse_ms && fUse_eionisation) {
486  n_order++;
487  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
488  }
489  }
490 
491 
492  if(particleName == "adj_gamma") {
493  G4int n_order=0;
494  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
496  pmanager->AddProcess(theAlongStepWeightCorrection);
497  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
498  idxAlongStep,1);
499 
500  if (fUse_brem && fUse_eionisation) {
501  pmanager->AddProcess(theeInverseBremsstrahlungProdToProjCase);
502  n_order++;
503  pmanager->SetProcessOrdering(
504  theeInverseBremsstrahlungProdToProjCase,
505  idxPostStep,n_order);
506  }
507  if (fUse_compton) {
508  pmanager->AddDiscreteProcess(theeInverseComptonProjToProjCase);
509  n_order++;
510  pmanager->SetProcessOrdering(theeInverseComptonProjToProjCase,
511  idxPostStep,n_order);
512  }
513  }
514 
515  if (particleName == "gamma") {
516  if (fUse_compton) {
517  pmanager->AddDiscreteProcess(theComptonScattering);
519  RegisterEmProcess(theComptonScattering,particle);
520  }
521  if (fUse_peeffect) {
522  pmanager->AddDiscreteProcess(thePEEffect);
524  RegisterEmProcess(thePEEffect,particle);
525  }
526  if (fUse_gamma_conversion) {
527  pmanager->AddDiscreteProcess(theGammaConversion);
528  }
529  }
530 
531  if (particleName == "e+" && fUse_gamma_conversion) {//positron
532  G4VProcess* theeplusMultipleScattering = new G4eMultipleScattering();
533  G4VProcess* theeplusIonisation = new G4eIonisation();
534  G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
535  G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
536 
537  // add processes
538  pmanager->AddProcess(theeplusMultipleScattering);
539  pmanager->AddProcess(theeplusIonisation);
540  pmanager->AddProcess(theeplusBremsstrahlung);
541  pmanager->AddProcess(theeplusAnnihilation);
542 
543  // set ordering for AtRestDoIt
544  pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
545 
546  // set ordering for AlongStepDoIt
547  pmanager->SetProcessOrdering(theeplusMultipleScattering,
548  idxAlongStep,1);
549  pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep,2);
550  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxAlongStep,3);
551 
552  // set ordering for PostStepDoIt
553  pmanager->SetProcessOrdering(theeplusMultipleScattering,
554  idxPostStep,1);
555  pmanager->SetProcessOrdering(theeplusIonisation,idxPostStep,2);
556  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxPostStep,3);
557  pmanager->SetProcessOrdering(theeplusAnnihilation,idxPostStep,4);
558  }
559  if (particleName == "proton" && fUse_pionisation) {
560  if (fUse_ms && fUse_pionisation) pmanager->AddProcess(thepMS);
561 
562  if (fUse_pionisation){
563  pmanager->AddProcess(fPIonisation);
565  RegisterEnergyLossProcess(fPIonisation,particle);
566  }
567 
568  G4int n_order=0;
569  if (fUse_ms && fUse_pionisation) {
570  n_order++;
571  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
572  }
573 
574  if (fUse_pionisation) {
575  n_order++;
576  pmanager->SetProcessOrdering(fPIonisation,idxAlongStep,n_order);
577  }
578 
579  n_order=0;
580  if (fUse_ms && fUse_pionisation) {
581  n_order++;
582  pmanager->SetProcessOrdering(thepMS, idxPostStep,n_order);
583  }
584 
585  if (fUse_pionisation) {
586  n_order++;
587  pmanager->SetProcessOrdering(fPIonisation,idxPostStep,n_order);
588  }
589 
590  }
591 
592  if (particleName == "adj_proton" && fUse_pionisation) {
593  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
594  if (fUse_pionisation ) {
595  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
596  theContinuousGainOfEnergy->SetLossFluctuations(
598  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(fPIonisation);
599  theContinuousGainOfEnergy->SetDirectParticle(G4Proton::Proton());
600  pmanager->AddProcess(theContinuousGainOfEnergy);
601  }
602 
603  G4int n_order=0;
604  if (fUse_ms) {
605  n_order++;
606  pmanager->AddProcess(thepMS);
607  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
608  }
609 
610  n_order++;
611  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
612  n_order);
613 
614  n_order++;
615  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
617  pmanager->AddProcess(theAlongStepWeightCorrection);
618  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
619  idxAlongStep,
620  n_order);
621  n_order=0;
622  if (fUse_pionisation) {
623  pmanager->AddProcess(thepInverseIonisationProjToProjCase);
624  n_order++;
625  pmanager->SetProcessOrdering(
626  thepInverseIonisationProjToProjCase,
627  idxPostStep,n_order);
628  }
629 
630  if (fUse_ms && fUse_pionisation) {
631  n_order++;
632  pmanager->SetProcessOrdering(thepMS,idxPostStep,n_order);
633  }
634  }
635  }
636 }
637 
638 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
639 
640 #include "G4Decay.hh"
642 {
643  // Add Decay Process
644  G4Decay* theDecayProcess = new G4Decay();
646  particleIterator->reset();
647  while( (*particleIterator)() ){
648  G4ParticleDefinition* particle = particleIterator->value();
649  G4ProcessManager* pmanager = particle->GetProcessManager();
650  if (theDecayProcess->IsApplicable(*particle)) {
651  pmanager ->AddProcess(theDecayProcess);
652  // set ordering for PostStepDoIt and AtRestDoIt
653  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
654  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
655  }
656  }
657 }
658 
659 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
660 
662 {
663  if (verboseLevel >0){
664  G4cout << "G4AdjointPhysicsList::SetCuts:";
665  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
666  }
667 
668  // set cut values for gamma at first and for e- second and next for e+,
669  // because some processes for e+/e- need cut values for gamma
670  //
671  SetCutValue(defaultCutValue, "gamma");
674 
676 }
677 
678 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
static G4AdjointSimManager * GetInstance()
static G4AdjointGamma * AdjointGamma()
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static const double MeV
Definition: G4SIunits.hh:211
static G4AdjointGamma * AdjointGammaDefinition()
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
G4AdjointPhysicsMessenger * fPhysicsMessenger
void SetProcessManager(G4ProcessManager *aProcessManager)
void SetCutValue(G4double aCut, const G4String &pname)
void SetProcessOrderingToFirst(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
virtual void ConstructParticle()
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4KaonZero * KaonZeroDefinition()
Definition: G4KaonZero.cc:99
static G4AntiKaonZero * AntiKaonZeroDefinition()
static G4AdjointElectron * AdjointElectron()
static G4KaonZeroShort * KaonZeroShortDefinition()
G4ProcessManager * GetProcessManager() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
void SetLowEnergyLimit(G4double aVal)
int G4int
Definition: G4Types.hh:78
static G4AntiNeutron * AntiNeutronDefinition()
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:103
static G4AdjointProton * AdjointProtonDefinition()
void SetHighEnergyLimit(G4double aVal)
void RegisterAdjointParticle(G4ParticleDefinition *aPartDef)
void SetLossFluctuationFlag(bool aBool)
void SetDirectEnergyLossProcess(G4VEnergyLossProcess *aProcess)
const G4String & GetParticleName() const
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
void SetLossFluctuations(G4bool val)
static G4KaonZeroLong * KaonZeroLongDefinition()
void SetVerboseLevel(G4int value)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
static G4NeutrinoE * NeutrinoEDefinition()
Definition: G4NeutrinoE.cc:80
void SetDirectParticle(G4ParticleDefinition *p)
static G4AntiNeutrinoMu * AntiNeutrinoMuDefinition()
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetUseMatrix(G4bool aBool)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
void SetDirectProcess(G4VEmProcess *aProcess)
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4EtaPrime * EtaPrimeDefinition()
Definition: G4EtaPrime.cc:106
static G4AdjointProton * AdjointProton()
Definition of the G4AdjointPhysicsList class.
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
static G4AdjointElectron * AdjointElectronDefinition()
void ConsiderParticleAsPrimary(const G4String &particle_name)
virtual void SetCSBiasingFactor(G4double aVal)
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
G4eIonisation * fEminusIonisation
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
static G4ChargedGeantino * ChargedGeantinoDefinition()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
static G4AntiNeutrinoE * AntiNeutrinoEDefinition()
static G4OpticalPhoton * OpticalPhotonDefinition()
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80
static const double mm
Definition: G4SIunits.hh:114
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4AdjointCSManager * GetAdjointCSManager()
Definition of the G4AdjointPhysicsMessenger class.
G4hIonisation * fPIonisation
static G4Eta * EtaDefinition()
Definition: G4Eta.cc:104
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81