Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/biasing/ReverseMC01/src/G4AdjointPhysicsList.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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