ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4AdjointPhysicsList.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4AdjointPhysicsList.cc
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 //
31 // Class Name: G4AdjointPhysicsList
32 // Author: L. Desorgher
33 // Organisation: SpaceIT GmbH
34 // Contract: ESA contract 21435/08/NL/AT
35 // Customer: ESA/ESTEC
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40 
41 #include "G4AdjointPhysicsList.hh"
42 #include "G4ProcessManager.hh"
43 #include "G4ParticleTypes.hh"
45 #include "G4SystemOfUnits.hh"
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48 
51  fEminusIonisation(0),fPIonisation(0),
52  fUse_forced_interaction(true),
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"
184 #include "G4hMultipleScattering.hh"
185 #include "G4eIonisation.hh"
186 #include "G4eBremsstrahlung.hh"
187 #include "G4eplusAnnihilation.hh"
188 #include "G4hIonisation.hh"
189 #include "G4ionIonisation.hh"
190 //#include "G4IonParametrisedLossModel.hh"
191 
192 #include "G4eBremsstrahlung.hh"
194 #include "G4eInverseIonisation.hh"
196 #include "G4AdjointCSManager.hh"
199 #include "G4AdjointComptonModel.hh"
200 #include "G4eInverseCompton.hh"
201 #include "G4InversePEEffect.hh"
204 #include "G4hInverseIonisation.hh"
207 #include "G4IonInverseIonisation.hh"
209 
210 #include "G4AdjointSimManager.hh"
212 
213 #include "G4SystemOfUnits.hh"
214 #include "G4PhysicalConstants.hh"
215 #include "G4UrbanMscModel.hh"
216 #include "G4UrbanAdjointMscModel.hh"
217 #include "G4UrbanAdjointMscModel.hh"
219 
220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
221 
223 {
224  G4AdjointCSManager* theCSManager =
226  G4AdjointSimManager* theAdjointSimManager =
228 
229  theCSManager->RegisterAdjointParticle(
231 
233  theCSManager->RegisterAdjointParticle(
235 
236  if (fUse_eionisation) {
238  new G4eIonisation();
240  }
241  if (fUse_pionisation) {
244  theCSManager->RegisterAdjointParticle(
246  }
247 
248  G4eBremsstrahlung* theeminusBremsstrahlung = 0;
250  theeminusBremsstrahlung = new G4eBremsstrahlung();
251 
252  G4ComptonScattering* theComptonScattering =0;
253  if (fUse_compton) theComptonScattering = new G4ComptonScattering();
254 
255  G4PhotoElectricEffect* thePEEffect =0;
256  if (fUse_peeffect) thePEEffect = new G4PhotoElectricEffect();
257 
258  G4eMultipleScattering* theeminusMS = 0;
259  G4hMultipleScattering* thepMS= 0;
260  G4eAdjointMultipleScattering* theeminusAdjointMS = 0;
261  if (fUse_ms) {
262  theeminusMS = new G4eMultipleScattering();
263  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
264  theeminusMS->SetEmModel(msc1);
265  theeminusAdjointMS = new G4eAdjointMultipleScattering();
267  theeminusAdjointMS->SetEmModel(msc2);
268  thepMS = new G4hMultipleScattering();
269  }
270 
271  G4VProcess* theGammaConversion =0;
272  if (fUse_gamma_conversion) theGammaConversion = new G4GammaConversion();
273  //Define adjoint e- ionisation
274  //-------------------
275  G4AdjointeIonisationModel* theeInverseIonisationModel = 0;
276  G4eInverseIonisation* theeInverseIonisationProjToProjCase = 0 ;
277  G4eInverseIonisation* theeInverseIonisationProdToProjCase = 0;
278  if (fUse_eionisation) {
279  theeInverseIonisationModel = new G4AdjointeIonisationModel();
280  theeInverseIonisationModel->SetHighEnergyLimit(
282  theeInverseIonisationModel->SetLowEnergyLimit(
284  theeInverseIonisationModel->SetCSBiasingFactor(
286  theeInverseIonisationProjToProjCase =
287  new G4eInverseIonisation(true,"Inv_eIon",
288  theeInverseIonisationModel);
289  theeInverseIonisationProdToProjCase =
290  new G4eInverseIonisation(false,"Inv_eIon1",
291  theeInverseIonisationModel);
292  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
293  }
294 
295  //Define adjoint Bremsstrahlung
296  //-------------------------------
297  G4AdjointBremsstrahlungModel* theeInverseBremsstrahlungModel = 0;
298  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProjToProjCase = 0;
299  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProdToProjCase = 0;
300  G4AdjointForcedInteractionForGamma* theForcedInteractionForGamma = 0;
301  if (fUse_brem && fUse_eionisation) {
302  theeInverseBremsstrahlungModel = new G4AdjointBremsstrahlungModel();
303  theeInverseBremsstrahlungModel->SetHighEnergyLimit(fEmax_adj_models*1.01);
304  theeInverseBremsstrahlungModel->SetLowEnergyLimit(fEmin_adj_models);
305  theeInverseBremsstrahlungModel->SetCSBiasingFactor(
307  theeInverseBremsstrahlungProjToProjCase = new G4eInverseBremsstrahlung(
308  true,"Inv_eBrem",theeInverseBremsstrahlungModel);
309  theeInverseBremsstrahlungProdToProjCase = new G4eInverseBremsstrahlung(
310  false,"Inv_eBrem1",theeInverseBremsstrahlungModel);
311  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
312  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
313 
314  if (!fUse_forced_interaction) theeInverseBremsstrahlungProdToProjCase
315  = new G4eInverseBremsstrahlung(false,G4String("Inv_eBrem1"),
316  theeInverseBremsstrahlungModel);
317  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
318  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
320  theForcedInteractionForGamma =
321  new G4AdjointForcedInteractionForGamma("ReverseGammaForcedInteraction");
322  theForcedInteractionForGamma->RegisterAdjointBremModel(
323  theeInverseBremsstrahlungModel);
324  }
325  }
326 
327 
328  //Define adjoint Compton
329  //---------------------
330 
331  G4AdjointComptonModel* theeInverseComptonModel = 0;
332  G4eInverseCompton* theeInverseComptonProjToProjCase = 0;
333  G4eInverseCompton* theeInverseComptonProdToProjCase = 0;
334 
335  if (fUse_compton) {
336  theeInverseComptonModel = new G4AdjointComptonModel();
337  theeInverseComptonModel->SetHighEnergyLimit(fEmax_adj_models);
338  theeInverseComptonModel->SetLowEnergyLimit(fEmin_adj_models);
339  theeInverseComptonModel->SetDirectProcess(theComptonScattering);
340  theeInverseComptonModel->SetUseMatrix(false);
341 
342  theeInverseComptonModel->SetCSBiasingFactor( fCS_biasing_factor_compton);
343  if (!fUse_forced_interaction) theeInverseComptonProjToProjCase =
344  new G4eInverseCompton(true,"Inv_Compt",theeInverseComptonModel);
345  theeInverseComptonProdToProjCase = new G4eInverseCompton(false,"Inv_Compt1",
346  theeInverseComptonModel);
348  if (!theForcedInteractionForGamma ) theForcedInteractionForGamma =
349  new G4AdjointForcedInteractionForGamma("ReverseGammaForcedInteraction");
350  theForcedInteractionForGamma->
351  RegisterAdjointComptonModel(theeInverseComptonModel);
352  }
353  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
354  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
355  }
356 
357  //Define adjoint PEEffect
358  //---------------------
359  G4AdjointPhotoElectricModel* theInversePhotoElectricModel = 0;
360  G4InversePEEffect* theInversePhotoElectricProcess = 0;
361 
362  if (fUse_peeffect) {
363  theInversePhotoElectricModel = new G4AdjointPhotoElectricModel();
364  theInversePhotoElectricModel->SetHighEnergyLimit(fEmax_adj_models);
365  theInversePhotoElectricModel->SetLowEnergyLimit(fEmin_adj_models);
366  theInversePhotoElectricModel->SetCSBiasingFactor(
368  theInversePhotoElectricProcess = new G4InversePEEffect("Inv_PEEffect",
369  theInversePhotoElectricModel);
370  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
371  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
372  }
373 
374 
375  //Define adjoint ionisation for protons
376  //---------------------
377  G4AdjointhIonisationModel* thepInverseIonisationModel = 0;
378  G4hInverseIonisation* thepInverseIonisationProjToProjCase = 0 ;
379  G4hInverseIonisation* thepInverseIonisationProdToProjCase = 0;
380  if (fUse_pionisation) {
381  thepInverseIonisationModel = new G4AdjointhIonisationModel(
382  G4Proton::Proton());
383  thepInverseIonisationModel->SetHighEnergyLimit(fEmax_adj_models);
384  thepInverseIonisationModel->SetLowEnergyLimit(fEmin_adj_models);
385  thepInverseIonisationModel->SetUseMatrix(false);
386  thepInverseIonisationProjToProjCase = new G4hInverseIonisation(true,
387  "Inv_pIon",thepInverseIonisationModel);
388  thepInverseIonisationProdToProjCase = new G4hInverseIonisation(false,
389  "Inv_pIon1",thepInverseIonisationModel);
390  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
391  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("proton"));
392  }
393 
394  //Declare the processes active for the different particles
395  //--------------------------------------------------------
397  particleIterator->reset();
398  while( (*particleIterator)() ){
400  G4ProcessManager* pmanager = particle->GetProcessManager();
401  if (!pmanager) {
402  pmanager = new G4ProcessManager(particle);
403  particle->SetProcessManager(pmanager);
404  }
405 
406  G4String particleName = particle->GetParticleName();
407  if (particleName == "e-") {
408  if (fUse_ms && fUse_eionisation) pmanager->AddProcess(theeminusMS);
409  if (fUse_eionisation){
410  pmanager->AddProcess(fEminusIonisation);
412  RegisterEnergyLossProcess(fEminusIonisation,particle);
413  }
414  if (fUse_brem && fUse_eionisation) {
415  pmanager->AddProcess(theeminusBremsstrahlung);
417  RegisterEnergyLossProcess(theeminusBremsstrahlung,particle);
418  }
419  G4int n_order=0;
420  if (fUse_ms && fUse_eionisation) {
421  n_order++;
422  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
423  }
424  if (fUse_eionisation) {
425  n_order++;
427  }
428  if (fUse_brem && fUse_eionisation) {
429  n_order++;
430  pmanager->SetProcessOrdering(theeminusBremsstrahlung,
431  idxAlongStep,n_order);
432  }
433  n_order=0;
434  if (fUse_ms && fUse_eionisation) {
435  n_order++;
436  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
437  }
438  if (fUse_eionisation) {
439  n_order++;
441  }
442  if (fUse_brem && fUse_eionisation) {
443  n_order++;
444  pmanager->SetProcessOrdering(theeminusBremsstrahlung,idxPostStep,
445  n_order);
446  }
447  }
448 
449  if (particleName == "adj_e-") {
450  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
451  if (fUse_eionisation ) {
452  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
453  theContinuousGainOfEnergy->SetLossFluctuations(
455  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(
457  theContinuousGainOfEnergy->SetDirectParticle(G4Electron::Electron());
458  pmanager->AddProcess(theContinuousGainOfEnergy);
459  }
460  G4int n_order=0;
461  if (fUse_ms) {
462  n_order++;
463  pmanager->AddProcess(theeminusAdjointMS);
464  pmanager->SetProcessOrdering(theeminusAdjointMS,
465  idxAlongStep,n_order);
466  }
467  n_order++;
468  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
469  n_order);
470 
471  n_order++;
472  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
474  pmanager->AddProcess(theAlongStepWeightCorrection);
475  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
476  idxAlongStep,
477  n_order);
478  n_order=0;
479  if (fUse_eionisation) {
480  pmanager->AddProcess(theeInverseIonisationProjToProjCase);
481  pmanager->AddProcess(theeInverseIonisationProdToProjCase);
482  n_order++;
483  pmanager->SetProcessOrdering(theeInverseIonisationProjToProjCase,
484  idxPostStep,n_order);
485  n_order++;
486  pmanager->SetProcessOrdering(theeInverseIonisationProdToProjCase,
487  idxPostStep,n_order);
488  }
489  if (fUse_brem && fUse_eionisation) {
490  pmanager->AddProcess(theeInverseBremsstrahlungProjToProjCase);
491  n_order++;
492  pmanager->SetProcessOrdering(
493  theeInverseBremsstrahlungProjToProjCase,
494  idxPostStep,n_order);
495  }
496 
497  if (fUse_compton) {
498  pmanager->AddProcess(theeInverseComptonProdToProjCase);
499  n_order++;
500  pmanager->SetProcessOrdering(theeInverseComptonProdToProjCase,
501  idxPostStep,n_order);
502  }
503  if (fUse_peeffect) {
504  pmanager->AddDiscreteProcess(theInversePhotoElectricProcess);
505  n_order++;
506  pmanager->SetProcessOrdering(theInversePhotoElectricProcess,
507  idxPostStep,n_order);
508  }
509  if (fUse_pionisation) {
510  pmanager->AddProcess(thepInverseIonisationProdToProjCase);
511  n_order++;
512  pmanager->SetProcessOrdering(thepInverseIonisationProdToProjCase,
513  idxPostStep,n_order);
514  }
515  if (fUse_ms && fUse_eionisation) {
516  n_order++;
517  pmanager->SetProcessOrdering(theeminusAdjointMS,
518  idxPostStep,n_order);
519  }
520  }
521 
522 
523  if(particleName == "adj_gamma") {
524  G4int n_order=0;
526  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
528  pmanager->AddProcess(theAlongStepWeightCorrection);
529  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
530  idxAlongStep,1);
531 
532  if (fUse_brem && fUse_eionisation) {
533  pmanager->AddProcess(theeInverseBremsstrahlungProdToProjCase);
534  n_order++;
535  pmanager->SetProcessOrdering(
536  theeInverseBremsstrahlungProdToProjCase,
537  idxPostStep,n_order);
538  }
539  if (fUse_compton) {
540  pmanager->AddDiscreteProcess(theeInverseComptonProjToProjCase);
541  n_order++;
542  pmanager->SetProcessOrdering(theeInverseComptonProjToProjCase,
543  idxPostStep,n_order);
544  }
545  }
546  else {
547  if (theForcedInteractionForGamma) {
548  pmanager->AddProcess(theForcedInteractionForGamma);
549  n_order++;
550  pmanager->SetProcessOrdering(theForcedInteractionForGamma,
551  idxPostStep,n_order);
552  pmanager->SetProcessOrdering(theForcedInteractionForGamma,
553  idxAlongStep,n_order);
554  }
555  }
556  }
557 
558  if (particleName == "gamma") {
559  if (fUse_compton) {
560  pmanager->AddDiscreteProcess(theComptonScattering);
562  RegisterEmProcess(theComptonScattering,particle);
563  }
564  if (fUse_peeffect) {
565  pmanager->AddDiscreteProcess(thePEEffect);
567  RegisterEmProcess(thePEEffect,particle);
568  }
569  if (fUse_gamma_conversion) {
570  pmanager->AddDiscreteProcess(theGammaConversion);
571  }
572  }
573 
574  if (particleName == "e+" && fUse_gamma_conversion) {//positron
575  G4VProcess* theeplusMultipleScattering = new G4eMultipleScattering();
576  G4VProcess* theeplusIonisation = new G4eIonisation();
577  G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
578  G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
579 
580  // add processes
581  pmanager->AddProcess(theeplusMultipleScattering);
582  pmanager->AddProcess(theeplusIonisation);
583  pmanager->AddProcess(theeplusBremsstrahlung);
584  pmanager->AddProcess(theeplusAnnihilation);
585 
586  // set ordering for AtRestDoIt
587  pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
588 
589  // set ordering for AlongStepDoIt
590  pmanager->SetProcessOrdering(theeplusMultipleScattering,
591  idxAlongStep,1);
592  pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep,2);
593  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxAlongStep,3);
594 
595  // set ordering for PostStepDoIt
596  pmanager->SetProcessOrdering(theeplusMultipleScattering,
597  idxPostStep,1);
598  pmanager->SetProcessOrdering(theeplusIonisation,idxPostStep,2);
599  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxPostStep,3);
600  pmanager->SetProcessOrdering(theeplusAnnihilation,idxPostStep,4);
601  }
602  if (particleName == "proton" && fUse_pionisation) {
603  if (fUse_ms && fUse_pionisation) pmanager->AddProcess(thepMS);
604 
605  if (fUse_pionisation){
606  pmanager->AddProcess(fPIonisation);
608  RegisterEnergyLossProcess(fPIonisation,particle);
609  }
610 
611  G4int n_order=0;
612  if (fUse_ms && fUse_pionisation) {
613  n_order++;
614  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
615  }
616 
617  if (fUse_pionisation) {
618  n_order++;
619  pmanager->SetProcessOrdering(fPIonisation,idxAlongStep,n_order);
620  }
621 
622  n_order=0;
623  if (fUse_ms && fUse_pionisation) {
624  n_order++;
625  pmanager->SetProcessOrdering(thepMS, idxPostStep,n_order);
626  }
627 
628  if (fUse_pionisation) {
629  n_order++;
630  pmanager->SetProcessOrdering(fPIonisation,idxPostStep,n_order);
631  }
632 
633  }
634 
635  if (particleName == "adj_proton" && fUse_pionisation) {
636  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
637  if (fUse_pionisation ) {
638  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
639  theContinuousGainOfEnergy->SetLossFluctuations(
641  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(fPIonisation);
642  theContinuousGainOfEnergy->SetDirectParticle(G4Proton::Proton());
643  pmanager->AddProcess(theContinuousGainOfEnergy);
644  }
645 
646  G4int n_order=0;
647  if (fUse_ms) {
648  n_order++;
649  pmanager->AddProcess(thepMS);
650  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
651  }
652 
653  n_order++;
654  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
655  n_order);
656 
657  n_order++;
658  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
660  pmanager->AddProcess(theAlongStepWeightCorrection);
661  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
662  idxAlongStep,
663  n_order);
664  n_order=0;
665  if (fUse_pionisation) {
666  pmanager->AddProcess(thepInverseIonisationProjToProjCase);
667  n_order++;
668  pmanager->SetProcessOrdering(
669  thepInverseIonisationProjToProjCase,
670  idxPostStep,n_order);
671  }
672 
673  if (fUse_ms && fUse_pionisation) {
674  n_order++;
675  pmanager->SetProcessOrdering(thepMS,idxPostStep,n_order);
676  }
677  }
678  }
679 }
680 
681 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
682 
683 #include "G4Decay.hh"
685 {
686  // Add Decay Process
687  G4Decay* theDecayProcess = new G4Decay();
689  particleIterator->reset();
690  while( (*particleIterator)() ){
692  G4ProcessManager* pmanager = particle->GetProcessManager();
693  if (theDecayProcess->IsApplicable(*particle)) {
694  pmanager ->AddProcess(theDecayProcess);
695  // set ordering for PostStepDoIt and AtRestDoIt
696  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
697  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
698  }
699  }
700 }
701 
702 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
703 
705 {
706  if (verboseLevel >0){
707  G4cout << "G4AdjointPhysicsList::SetCuts:";
708  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
709  }
710 
711  // set cut values for gamma at first and for e- second and next for e+,
712  // because some processes for e+/e- need cut values for gamma
713  //
714  SetCutValue(defaultCutValue, "gamma");
717 
719 }
720 
721 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......