Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
26 /// \file biasing/ReverseMC01/src/G4AdjointPhysicsList.cc
27 /// \brief Implementation of the G4AdjointPhysicsList class
28 //
29 // $Id: G4AdjointPhysicsList.cc 70966 2013-06-07 15:22:09Z gcosmo $
30 //
31 //////////////////////////////////////////////////////////////
32 // Class Name: G4AdjointPhysicsList
33 // Author: L. Desorgher
34 // Organisation: SpaceIT GmbH
35 // Contract: ESA contract 21435/08/NL/AT
36 // Customer: ESA/ESTEC
37 //////////////////////////////////////////////////////////////
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);
62  fPhysicsMessenger = new G4AdjointPhysicsMessenger(this);
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 
221  if (fUse_brem || fUse_peeffect ||fUse_compton)
222  theCSManager->RegisterAdjointParticle(
224 
225  if (fUse_eionisation) {
227  new G4eIonisation();
228  fEminusIonisation->SetLossFluctuations(fUse_egain_fluctuation);
229  }
230 
231  if (fUse_pionisation) {
233  fPIonisation->SetLossFluctuations(fUse_egain_fluctuation);
234  theCSManager->RegisterAdjointParticle(
236  }
237 
238  G4eBremsstrahlung* theeminusBremsstrahlung = 0;
239  if (fUse_brem && fUse_eionisation)
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(
266  fEmax_adj_models);
267  theeInverseIonisationModel->SetLowEnergyLimit(
268  fEmin_adj_models);
269  theeInverseIonisationModel->SetCSBiasingFactor(
270  fCS_biasing_factor_ionisation);
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(
291  fCS_biasing_factor_brem);
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(
333  fCS_biasing_factor_PEeffect);
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  //--------------------------------------------------------
364  theParticleIterator->reset();
365  while( (*theParticleIterator)() ){
366  G4ParticleDefinition* particle = theParticleIterator->value();
367  G4ProcessManager* pmanager = particle->GetProcessManager();
368  if (!pmanager) {
369  pmanager = new G4ProcessManager(particle);
370  particle->SetProcessManager(pmanager);
371  }
372 
373  G4String particleName = particle->GetParticleName();
374  if (particleName == "e-") {
375  if (fUse_ms && fUse_eionisation) pmanager->AddProcess(theeminusMS);
376  if (fUse_eionisation){
377  pmanager->AddProcess(fEminusIonisation);
379  RegisterEnergyLossProcess(fEminusIonisation,particle);
380  }
381  if (fUse_brem && fUse_eionisation) {
382  pmanager->AddProcess(theeminusBremsstrahlung);
384  RegisterEnergyLossProcess(theeminusBremsstrahlung,particle);
385  }
386  G4int n_order=0;
387  if (fUse_ms && fUse_eionisation) {
388  n_order++;
389  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
390  }
391  if (fUse_eionisation) {
392  n_order++;
394  }
395  if (fUse_brem && fUse_eionisation) {
396  n_order++;
397  pmanager->SetProcessOrdering(theeminusBremsstrahlung,
398  idxAlongStep,n_order);
399  }
400  n_order=0;
401  if (fUse_ms && fUse_eionisation) {
402  n_order++;
403  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
404  }
405  if (fUse_eionisation) {
406  n_order++;
408  }
409  if (fUse_brem && fUse_eionisation) {
410  n_order++;
411  pmanager->SetProcessOrdering(theeminusBremsstrahlung,idxPostStep,
412  n_order);
413  }
414  }
415 
416  if (particleName == "adj_e-") {
417  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
418  if (fUse_eionisation ) {
419  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
420  theContinuousGainOfEnergy->SetLossFluctuations(
421  fUse_egain_fluctuation);
422  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(
424  theContinuousGainOfEnergy->SetDirectParticle(G4Electron::Electron());
425  pmanager->AddProcess(theContinuousGainOfEnergy);
426  }
427  G4int n_order=0;
428  if (fUse_ms) {
429  n_order++;
430  pmanager->AddProcess(theeminusMS);
431  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
432  }
433 
434  n_order++;
435  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
436  n_order);
437 
438 
439  n_order++;
440  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
442  pmanager->AddProcess(theAlongStepWeightCorrection);
443  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
444  idxAlongStep,
445  n_order);
446  n_order=0;
447  if (fUse_eionisation) {
448  pmanager->AddProcess(theeInverseIonisationProjToProjCase);
449  pmanager->AddProcess(theeInverseIonisationProdToProjCase);
450  n_order++;
451  pmanager->SetProcessOrdering(theeInverseIonisationProjToProjCase,
452  idxPostStep,n_order);
453  n_order++;
454  pmanager->SetProcessOrdering(theeInverseIonisationProdToProjCase,
455  idxPostStep,n_order);
456  }
457 
458  if (fUse_brem && fUse_eionisation) {
459  pmanager->AddProcess(theeInverseBremsstrahlungProjToProjCase);
460  n_order++;
461  pmanager->SetProcessOrdering(
462  theeInverseBremsstrahlungProjToProjCase,
463  idxPostStep,n_order);
464  }
465 
466  if (fUse_compton) {
467  pmanager->AddProcess(theeInverseComptonProdToProjCase);
468  n_order++;
469  pmanager->SetProcessOrdering(theeInverseComptonProdToProjCase,
470  idxPostStep,n_order);
471  }
472  if (fUse_peeffect) {
473  pmanager->AddDiscreteProcess(theInversePhotoElectricProcess);
474  n_order++;
475  pmanager->SetProcessOrdering(theInversePhotoElectricProcess,
476  idxPostStep,n_order);
477  }
478  if (fUse_pionisation) {
479  pmanager->AddProcess(thepInverseIonisationProdToProjCase);
480  n_order++;
481  pmanager->SetProcessOrdering(thepInverseIonisationProdToProjCase,
482  idxPostStep,n_order);
483  }
484  if (fUse_ms && fUse_eionisation) {
485  n_order++;
486  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
487  }
488  }
489 
490 
491  if(particleName == "adj_gamma") {
492  G4int n_order=0;
493  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
495  pmanager->AddProcess(theAlongStepWeightCorrection);
496  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
497  idxAlongStep,1);
498 
499  if (fUse_brem && fUse_eionisation) {
500  pmanager->AddProcess(theeInverseBremsstrahlungProdToProjCase);
501  n_order++;
502  pmanager->SetProcessOrdering(
503  theeInverseBremsstrahlungProdToProjCase,
504  idxPostStep,n_order);
505  }
506  if (fUse_compton) {
507  pmanager->AddDiscreteProcess(theeInverseComptonProjToProjCase);
508  n_order++;
509  pmanager->SetProcessOrdering(theeInverseComptonProjToProjCase,
510  idxPostStep,n_order);
511  }
512  }
513 
514  if (particleName == "gamma") {
515  if (fUse_compton) {
516  pmanager->AddDiscreteProcess(theComptonScattering);
518  RegisterEmProcess(theComptonScattering,particle);
519  }
520  if (fUse_peeffect) {
521  pmanager->AddDiscreteProcess(thePEEffect);
523  RegisterEmProcess(thePEEffect,particle);
524  }
525  if (fUse_gamma_conversion) {
526  pmanager->AddDiscreteProcess(theGammaConversion);
527  }
528  }
529 
530  if (particleName == "e+" && fUse_gamma_conversion) {//positron
531  G4VProcess* theeplusMultipleScattering = new G4eMultipleScattering();
532  G4VProcess* theeplusIonisation = new G4eIonisation();
533  G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
534  G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
535 
536  // add processes
537  pmanager->AddProcess(theeplusMultipleScattering);
538  pmanager->AddProcess(theeplusIonisation);
539  pmanager->AddProcess(theeplusBremsstrahlung);
540  pmanager->AddProcess(theeplusAnnihilation);
541 
542  // set ordering for AtRestDoIt
543  pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
544 
545  // set ordering for AlongStepDoIt
546  pmanager->SetProcessOrdering(theeplusMultipleScattering,
547  idxAlongStep,1);
548  pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep,2);
549  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxAlongStep,3);
550 
551  // set ordering for PostStepDoIt
552  pmanager->SetProcessOrdering(theeplusMultipleScattering,
553  idxPostStep,1);
554  pmanager->SetProcessOrdering(theeplusIonisation,idxPostStep,2);
555  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxPostStep,3);
556  pmanager->SetProcessOrdering(theeplusAnnihilation,idxPostStep,4);
557  }
558  if (particleName == "proton" && fUse_pionisation) {
559  if (fUse_ms && fUse_pionisation) pmanager->AddProcess(thepMS);
560 
561  if (fUse_pionisation){
562  pmanager->AddProcess(fPIonisation);
564  RegisterEnergyLossProcess(fPIonisation,particle);
565  }
566 
567  G4int n_order=0;
568  if (fUse_ms && fUse_pionisation) {
569  n_order++;
570  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
571  }
572 
573  if (fUse_pionisation) {
574  n_order++;
575  pmanager->SetProcessOrdering(fPIonisation,idxAlongStep,n_order);
576  }
577 
578  n_order=0;
579  if (fUse_ms && fUse_pionisation) {
580  n_order++;
581  pmanager->SetProcessOrdering(thepMS, idxPostStep,n_order);
582  }
583 
584  if (fUse_pionisation) {
585  n_order++;
586  pmanager->SetProcessOrdering(fPIonisation,idxPostStep,n_order);
587  }
588 
589  }
590 
591  if (particleName == "adj_proton" && fUse_pionisation) {
592  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
593  if (fUse_pionisation ) {
594  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
595  theContinuousGainOfEnergy->SetLossFluctuations(
596  fUse_egain_fluctuation);
597  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(fPIonisation);
598  theContinuousGainOfEnergy->SetDirectParticle(G4Proton::Proton());
599  pmanager->AddProcess(theContinuousGainOfEnergy);
600  }
601 
602  G4int n_order=0;
603  if (fUse_ms) {
604  n_order++;
605  pmanager->AddProcess(thepMS);
606  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
607  }
608 
609  n_order++;
610  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
611  n_order);
612 
613  n_order++;
614  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
616  pmanager->AddProcess(theAlongStepWeightCorrection);
617  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
618  idxAlongStep,
619  n_order);
620  n_order=0;
621  if (fUse_pionisation) {
622  pmanager->AddProcess(thepInverseIonisationProjToProjCase);
623  n_order++;
624  pmanager->SetProcessOrdering(
625  thepInverseIonisationProjToProjCase,
626  idxPostStep,n_order);
627  }
628 
629  if (fUse_ms && fUse_pionisation) {
630  n_order++;
631  pmanager->SetProcessOrdering(thepMS,idxPostStep,n_order);
632  }
633  }
634  }
635 }
636 
637 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
638 
639 #include "G4Decay.hh"
641 {
642  // Add Decay Process
643  G4Decay* theDecayProcess = new G4Decay();
644  theParticleIterator->reset();
645  while( (*theParticleIterator)() ){
646  G4ParticleDefinition* particle = theParticleIterator->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 
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");
672 
674 }
675 
676 //....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 G4AdjointGamma * AdjointGammaDefinition()
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
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()
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
void SetLowEnergyLimit(G4double aVal)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
static G4AntiNeutron * AntiNeutronDefinition()
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:99
const G4String & GetParticleName() const
static G4AdjointProton * AdjointProtonDefinition()
void SetHighEnergyLimit(G4double aVal)
void RegisterAdjointParticle(G4ParticleDefinition *aPartDef)
void SetLossFluctuationFlag(bool aBool)
void SetDirectEnergyLossProcess(G4VEnergyLossProcess *aProcess)
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:100
static G4AdjointProton * AdjointProton()
Definition of the G4AdjointPhysicsList class.
static G4AdjointElectron * AdjointElectronDefinition()
void ConsiderParticleAsPrimary(const G4String &particle_name)
virtual void SetCSBiasingFactor(G4double aVal)
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
G4eIonisation * fEminusIonisation
static G4ChargedGeantino * ChargedGeantinoDefinition()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4AntiNeutrinoE * AntiNeutrinoEDefinition()
static G4OpticalPhoton * OpticalPhotonDefinition()
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4AdjointCSManager * GetAdjointCSManager()
#define theParticleIterator
Definition of the G4AdjointPhysicsMessenger class.
G4hIonisation * fPIonisation
static G4Eta * EtaDefinition()
Definition: G4Eta.cc:104
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81