Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4AdjointPhysicsList Class Reference

#include <G4AdjointPhysicsList.hh>

Inheritance diagram for G4AdjointPhysicsList:
Collaboration diagram for G4AdjointPhysicsList:

Public Member Functions

 G4AdjointPhysicsList ()
 
virtual ~G4AdjointPhysicsList ()
 
void SetLossFluctuationFlag (bool aBool)
 
void SetUseIonisation (bool aBool)
 
void SetUseProtonIonisation (bool aBool)
 
void SetUseBrem (bool aBool)
 
void SetUseCompton (bool aBool)
 
void SetUseMS (bool aBool)
 
void SetUsePEEffect (bool aBool)
 
void SetUseGammaConversion (bool aBool)
 
void SetUseEgainFluctuation (bool aBool)
 
void SetEminAdjModels (G4double aVal)
 
void SetEmaxAdjModels (G4double aVal)
 
- Public Member Functions inherited from G4VUserPhysicsList
 G4VUserPhysicsList ()
 
virtual ~G4VUserPhysicsList ()
 
 G4VUserPhysicsList (const G4VUserPhysicsList &)
 
G4VUserPhysicsListoperator= (const G4VUserPhysicsList &)
 
void Construct ()
 
void UseCoupledTransportation (G4bool vl=true)
 
void SetDefaultCutValue (G4double newCutValue)
 
G4double GetDefaultCutValue () const
 
void BuildPhysicsTable ()
 
void PreparePhysicsTable (G4ParticleDefinition *)
 
void BuildPhysicsTable (G4ParticleDefinition *)
 
G4bool StorePhysicsTable (const G4String &directory=".")
 
G4bool IsPhysicsTableRetrieved () const
 
G4bool IsStoredInAscii () const
 
const G4StringGetPhysicsTableDirectory () const
 
void SetPhysicsTableRetrieved (const G4String &directory="")
 
void SetStoredInAscii ()
 
void ResetPhysicsTableRetrieved ()
 
void ResetStoredInAscii ()
 
void DumpList () const
 
void DumpCutValuesTable (G4int flag=1)
 
void DumpCutValuesTableIfRequested ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void SetCutsWithDefault ()
 
void SetCutValue (G4double aCut, const G4String &pname)
 
G4double GetCutValue (const G4String &pname) const
 
void SetCutValue (G4double aCut, const G4String &pname, const G4String &rname)
 
void SetParticleCuts (G4double cut, G4ParticleDefinition *particle, G4Region *region=0)
 
void SetParticleCuts (G4double cut, const G4String &particleName, G4Region *region=0)
 
void SetCutsForRegion (G4double aCut, const G4String &rname)
 
void ResetCuts ()
 obsolete methods More...
 
void SetApplyCuts (G4bool value, const G4String &name)
 
G4bool GetApplyCuts (const G4String &name) const
 
void RemoveProcessManager ()
 
void AddProcessManager (G4ParticleDefinition *newParticle, G4ProcessManager *newManager=0)
 
void CheckParticleList ()
 
void DisableCheckParticleList ()
 
G4int GetInstanceID () const
 
void InitializeWorker ()
 

Protected Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
virtual void SetCuts ()
 
void ConstructBosons ()
 
void ConstructLeptons ()
 
void ConstructMesons ()
 
void ConstructBaryons ()
 
void ConstructAdjointParticles ()
 
void ConstructGeneral ()
 
void ConstructEM ()
 
- Protected Member Functions inherited from G4VUserPhysicsList
void AddTransportation ()
 
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
void BuildIntegralPhysicsTable (G4VProcess *, G4ParticleDefinition *)
 
virtual void RetrievePhysicsTable (G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
void InitializeProcessManager ()
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 

Protected Attributes

G4eIonisationfEminusIonisation
 
G4hIonisationfPIonisation
 
- Protected Attributes inherited from G4VUserPhysicsList
G4ParticleTabletheParticleTable
 
G4int verboseLevel
 
G4double defaultCutValue
 
G4bool isSetDefaultCutValue
 
G4ProductionCutsTablefCutsTable
 
G4bool fRetrievePhysicsTable
 
G4bool fStoredInAscii
 
G4bool fIsCheckedForRetrievePhysicsTable
 
G4bool fIsRestoredCutValues
 
G4String directoryPhysicsTable
 
G4bool fDisableCheckParticleList
 
G4int g4vuplInstanceID
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VUserPhysicsList
static const G4VUPLManagerGetSubInstanceManager ()
 
- Static Protected Attributes inherited from G4VUserPhysicsList
static G4RUN_DLL G4VUPLManager subInstanceManager
 

Detailed Description

Definition at line 58 of file G4AdjointPhysicsList.hh.

Constructor & Destructor Documentation

G4AdjointPhysicsList::G4AdjointPhysicsList ( )

Definition at line 50 of file G4AdjointPhysicsList.cc.

53  fUse_forced_interaction(true),
54  fUse_eionisation(true),fUse_pionisation(true),
55  fUse_brem(true),fUse_compton(true),fUse_ms(true),
56  fUse_egain_fluctuation(true),fUse_peeffect(true),
57  fEmin_adj_models(1.*keV), fEmax_adj_models(1.*MeV),
58  fCS_biasing_factor_compton(1.),fCS_biasing_factor_brem(1.),
59  fCS_biasing_factor_ionisation(1.),fCS_biasing_factor_PEeffect(1.)
60 {
61  defaultCutValue = 1.0*mm;
62  SetVerboseLevel(1);
63  fPhysicsMessenger = new G4AdjointPhysicsMessenger(this);
64 }
static constexpr double mm
Definition: G4SIunits.hh:115
void SetVerboseLevel(G4int value)
G4eIonisation * fEminusIonisation
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double keV
Definition: G4SIunits.hh:216
G4hIonisation * fPIonisation

Here is the call graph for this function:

G4AdjointPhysicsList::~G4AdjointPhysicsList ( )
virtual

Definition at line 68 of file G4AdjointPhysicsList.cc.

69 {
70 }

Member Function Documentation

void G4AdjointPhysicsList::ConstructAdjointParticles ( )
protected

Definition at line 156 of file G4AdjointPhysicsList.cc.

157 {
158 // adjoint_gammma
160 
161 // adjoint_electron
163 
164 // adjoint_proton
166 }
static G4AdjointGamma * AdjointGammaDefinition()
static G4AdjointProton * AdjointProtonDefinition()
static G4AdjointElectron * AdjointElectronDefinition()

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructBaryons ( )
protected

Definition at line 142 of file G4AdjointPhysicsList.cc.

143 {
144 // barions
149 }
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
static G4AntiNeutron * AntiNeutronDefinition()
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructBosons ( )
protected

Definition at line 93 of file G4AdjointPhysicsList.cc.

94 {
95  // pseudo-particles
98 
99  // gamma
101 
102  // optical photon
104 }
static G4Geantino * GeantinoDefinition()
Definition: G4Geantino.cc:82
static G4ChargedGeantino * ChargedGeantinoDefinition()
static G4OpticalPhoton * OpticalPhotonDefinition()
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructEM ( )
protected

Definition at line 223 of file G4AdjointPhysicsList.cc.

224 {
225  G4AdjointCSManager* theCSManager =
227  G4AdjointSimManager* theAdjointSimManager =
229 
230  theCSManager->RegisterAdjointParticle(
232 
233  if (fUse_brem || fUse_peeffect ||fUse_compton)
234  theCSManager->RegisterAdjointParticle(
236 
237  if (fUse_eionisation) {
239  new G4eIonisation();
240  fEminusIonisation->SetLossFluctuations(fUse_egain_fluctuation);
241  }
242  if (fUse_pionisation) {
244  fPIonisation->SetLossFluctuations(fUse_egain_fluctuation);
245  theCSManager->RegisterAdjointParticle(
247  }
248 
249  G4eBremsstrahlung* theeminusBremsstrahlung = 0;
250  if (fUse_brem && fUse_eionisation)
251  theeminusBremsstrahlung = new G4eBremsstrahlung();
252 
253  G4ComptonScattering* theComptonScattering =0;
254  if (fUse_compton) theComptonScattering = new G4ComptonScattering();
255 
256  G4PhotoElectricEffect* thePEEffect =0;
257  if (fUse_peeffect) thePEEffect = new G4PhotoElectricEffect();
258 
259  G4eMultipleScattering* theeminusMS = 0;
260  G4hMultipleScattering* thepMS= 0;
261  G4eAdjointMultipleScattering* theeminusAdjointMS = 0;
262  if (fUse_ms) {
263  theeminusMS = new G4eMultipleScattering();
264  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
265  msc1->SetNewDisplacementFlag(false);
266  theeminusMS->AddEmModel(0, msc1);
267  theeminusAdjointMS = new G4eAdjointMultipleScattering();
269  msc2->SetNewDisplacementFlag(false);
270  theeminusAdjointMS->AddEmModel(0, msc2);
271  thepMS = new G4hMultipleScattering();
272  }
273 
274  G4VProcess* theGammaConversion =0;
275  if (fUse_gamma_conversion) theGammaConversion = new G4GammaConversion();
276  //Define adjoint e- ionisation
277  //-------------------
278  G4AdjointeIonisationModel* theeInverseIonisationModel = 0;
279  G4eInverseIonisation* theeInverseIonisationProjToProjCase = 0 ;
280  G4eInverseIonisation* theeInverseIonisationProdToProjCase = 0;
281  if (fUse_eionisation) {
282  theeInverseIonisationModel = new G4AdjointeIonisationModel();
283  theeInverseIonisationModel->SetHighEnergyLimit(
284  fEmax_adj_models);
285  theeInverseIonisationModel->SetLowEnergyLimit(
286  fEmin_adj_models);
287  theeInverseIonisationModel->SetCSBiasingFactor(
288  fCS_biasing_factor_ionisation);
289  theeInverseIonisationProjToProjCase =
290  new G4eInverseIonisation(true,"Inv_eIon",
291  theeInverseIonisationModel);
292  theeInverseIonisationProdToProjCase =
293  new G4eInverseIonisation(false,"Inv_eIon1",
294  theeInverseIonisationModel);
295  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
296  }
297 
298  //Define adjoint Bremsstrahlung
299  //-------------------------------
300  G4AdjointBremsstrahlungModel* theeInverseBremsstrahlungModel = 0;
301  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProjToProjCase = 0;
302  G4eInverseBremsstrahlung* theeInverseBremsstrahlungProdToProjCase = 0;
303  G4AdjointForcedInteractionForGamma* theForcedInteractionForGamma = 0;
304  if (fUse_brem && fUse_eionisation) {
305  theeInverseBremsstrahlungModel = new G4AdjointBremsstrahlungModel();
306  theeInverseBremsstrahlungModel->SetHighEnergyLimit(fEmax_adj_models*1.01);
307  theeInverseBremsstrahlungModel->SetLowEnergyLimit(fEmin_adj_models);
308  theeInverseBremsstrahlungModel->SetCSBiasingFactor(
309  fCS_biasing_factor_brem);
310  theeInverseBremsstrahlungProjToProjCase = new G4eInverseBremsstrahlung(
311  true,"Inv_eBrem",theeInverseBremsstrahlungModel);
312  theeInverseBremsstrahlungProdToProjCase = new G4eInverseBremsstrahlung(
313  false,"Inv_eBrem1",theeInverseBremsstrahlungModel);
314  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
315  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
316 
317  if (!fUse_forced_interaction) theeInverseBremsstrahlungProdToProjCase
318  = new G4eInverseBremsstrahlung(false,G4String("Inv_eBrem1"),
319  theeInverseBremsstrahlungModel);
320  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
321  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
322  if (fUse_forced_interaction){
323  theForcedInteractionForGamma =
324  new G4AdjointForcedInteractionForGamma("ReverseGammaForcedInteraction");
325  theForcedInteractionForGamma->RegisterAdjointBremModel(
326  theeInverseBremsstrahlungModel);
327  }
328  }
329 
330 
331  //Define adjoint Compton
332  //---------------------
333 
334  G4AdjointComptonModel* theeInverseComptonModel = 0;
335  G4eInverseCompton* theeInverseComptonProjToProjCase = 0;
336  G4eInverseCompton* theeInverseComptonProdToProjCase = 0;
337 
338  if (fUse_compton) {
339  theeInverseComptonModel = new G4AdjointComptonModel();
340  theeInverseComptonModel->SetHighEnergyLimit(fEmax_adj_models);
341  theeInverseComptonModel->SetLowEnergyLimit(fEmin_adj_models);
342  theeInverseComptonModel->SetDirectProcess(theComptonScattering);
343  theeInverseComptonModel->SetUseMatrix(false);
344 
345  theeInverseComptonModel->SetCSBiasingFactor( fCS_biasing_factor_compton);
346  if (!fUse_forced_interaction) theeInverseComptonProjToProjCase =
347  new G4eInverseCompton(true,"Inv_Compt",theeInverseComptonModel);
348  theeInverseComptonProdToProjCase = new G4eInverseCompton(false,"Inv_Compt1",
349  theeInverseComptonModel);
350  if (fUse_forced_interaction){
351  if (!theForcedInteractionForGamma ) theForcedInteractionForGamma =
352  new G4AdjointForcedInteractionForGamma("ReverseGammaForcedInteraction");
353  theForcedInteractionForGamma->
354  RegisterAdjointComptonModel(theeInverseComptonModel);
355  }
356  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
357  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
358  }
359 
360  //Define adjoint PEEffect
361  //---------------------
362  G4AdjointPhotoElectricModel* theInversePhotoElectricModel = 0;
363  G4InversePEEffect* theInversePhotoElectricProcess = 0;
364 
365  if (fUse_peeffect) {
366  theInversePhotoElectricModel = new G4AdjointPhotoElectricModel();
367  theInversePhotoElectricModel->SetHighEnergyLimit(fEmax_adj_models);
368  theInversePhotoElectricModel->SetLowEnergyLimit(fEmin_adj_models);
369  theInversePhotoElectricModel->SetCSBiasingFactor(
370  fCS_biasing_factor_PEeffect);
371  theInversePhotoElectricProcess = new G4InversePEEffect("Inv_PEEffect",
372  theInversePhotoElectricModel);
373  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
374  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("gamma"));
375  }
376 
377 
378  //Define adjoint ionisation for protons
379  //---------------------
380  G4AdjointhIonisationModel* thepInverseIonisationModel = 0;
381  G4hInverseIonisation* thepInverseIonisationProjToProjCase = 0 ;
382  G4hInverseIonisation* thepInverseIonisationProdToProjCase = 0;
383  if (fUse_pionisation) {
384  thepInverseIonisationModel = new G4AdjointhIonisationModel(
385  G4Proton::Proton());
386  thepInverseIonisationModel->SetHighEnergyLimit(fEmax_adj_models);
387  thepInverseIonisationModel->SetLowEnergyLimit(fEmin_adj_models);
388  thepInverseIonisationModel->SetUseMatrix(false);
389  thepInverseIonisationProjToProjCase = new G4hInverseIonisation(true,
390  "Inv_pIon",thepInverseIonisationModel);
391  thepInverseIonisationProdToProjCase = new G4hInverseIonisation(false,
392  "Inv_pIon1",thepInverseIonisationModel);
393  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("e-"));
394  theAdjointSimManager->ConsiderParticleAsPrimary(G4String("proton"));
395  }
396 
397  //Declare the processes active for the different particles
398  //--------------------------------------------------------
400  particleIterator->reset();
401  while( (*particleIterator)() ){
402  G4ParticleDefinition* particle = particleIterator->value();
403  G4ProcessManager* pmanager = particle->GetProcessManager();
404  if (!pmanager) {
405  pmanager = new G4ProcessManager(particle);
406  particle->SetProcessManager(pmanager);
407  }
408 
409  G4String particleName = particle->GetParticleName();
410  if (particleName == "e-") {
411  if (fUse_ms && fUse_eionisation) pmanager->AddProcess(theeminusMS);
412  if (fUse_eionisation){
413  pmanager->AddProcess(fEminusIonisation);
415  RegisterEnergyLossProcess(fEminusIonisation,particle);
416  }
417  if (fUse_brem && fUse_eionisation) {
418  pmanager->AddProcess(theeminusBremsstrahlung);
420  RegisterEnergyLossProcess(theeminusBremsstrahlung,particle);
421  }
422  G4int n_order=0;
423  if (fUse_ms && fUse_eionisation) {
424  n_order++;
425  pmanager->SetProcessOrdering(theeminusMS, idxAlongStep,n_order);
426  }
427  if (fUse_eionisation) {
428  n_order++;
430  }
431  if (fUse_brem && fUse_eionisation) {
432  n_order++;
433  pmanager->SetProcessOrdering(theeminusBremsstrahlung,
434  idxAlongStep,n_order);
435  }
436  n_order=0;
437  if (fUse_ms && fUse_eionisation) {
438  n_order++;
439  pmanager->SetProcessOrdering(theeminusMS,idxPostStep,n_order);
440  }
441  if (fUse_eionisation) {
442  n_order++;
444  }
445  if (fUse_brem && fUse_eionisation) {
446  n_order++;
447  pmanager->SetProcessOrdering(theeminusBremsstrahlung,idxPostStep,
448  n_order);
449  }
450  }
451 
452  if (particleName == "adj_e-") {
453  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
454  if (fUse_eionisation ) {
455  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
456  theContinuousGainOfEnergy->SetLossFluctuations(
457  fUse_egain_fluctuation);
458  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(
460  theContinuousGainOfEnergy->SetDirectParticle(G4Electron::Electron());
461  pmanager->AddProcess(theContinuousGainOfEnergy);
462  }
463  G4int n_order=0;
464  if (fUse_ms) {
465  n_order++;
466  pmanager->AddProcess(theeminusAdjointMS);
467  pmanager->SetProcessOrdering(theeminusAdjointMS,
468  idxAlongStep,n_order);
469  }
470  n_order++;
471  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
472  n_order);
473 
474  n_order++;
475  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
477  pmanager->AddProcess(theAlongStepWeightCorrection);
478  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
479  idxAlongStep,
480  n_order);
481  n_order=0;
482  if (fUse_eionisation) {
483  pmanager->AddProcess(theeInverseIonisationProjToProjCase);
484  pmanager->AddProcess(theeInverseIonisationProdToProjCase);
485  n_order++;
486  pmanager->SetProcessOrdering(theeInverseIonisationProjToProjCase,
487  idxPostStep,n_order);
488  n_order++;
489  pmanager->SetProcessOrdering(theeInverseIonisationProdToProjCase,
490  idxPostStep,n_order);
491  }
492  if (fUse_brem && fUse_eionisation) {
493  pmanager->AddProcess(theeInverseBremsstrahlungProjToProjCase);
494  n_order++;
495  pmanager->SetProcessOrdering(
496  theeInverseBremsstrahlungProjToProjCase,
497  idxPostStep,n_order);
498  }
499 
500  if (fUse_compton) {
501  pmanager->AddProcess(theeInverseComptonProdToProjCase);
502  n_order++;
503  pmanager->SetProcessOrdering(theeInverseComptonProdToProjCase,
504  idxPostStep,n_order);
505  }
506  if (fUse_peeffect) {
507  pmanager->AddDiscreteProcess(theInversePhotoElectricProcess);
508  n_order++;
509  pmanager->SetProcessOrdering(theInversePhotoElectricProcess,
510  idxPostStep,n_order);
511  }
512  if (fUse_pionisation) {
513  pmanager->AddProcess(thepInverseIonisationProdToProjCase);
514  n_order++;
515  pmanager->SetProcessOrdering(thepInverseIonisationProdToProjCase,
516  idxPostStep,n_order);
517  }
518  if (fUse_ms && fUse_eionisation) {
519  n_order++;
520  pmanager->SetProcessOrdering(theeminusAdjointMS,
521  idxPostStep,n_order);
522  }
523  }
524 
525 
526  if(particleName == "adj_gamma") {
527  G4int n_order=0;
528  if (!fUse_forced_interaction){
529  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
531  pmanager->AddProcess(theAlongStepWeightCorrection);
532  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
533  idxAlongStep,1);
534 
535  if (fUse_brem && fUse_eionisation) {
536  pmanager->AddProcess(theeInverseBremsstrahlungProdToProjCase);
537  n_order++;
538  pmanager->SetProcessOrdering(
539  theeInverseBremsstrahlungProdToProjCase,
540  idxPostStep,n_order);
541  }
542  if (fUse_compton) {
543  pmanager->AddDiscreteProcess(theeInverseComptonProjToProjCase);
544  n_order++;
545  pmanager->SetProcessOrdering(theeInverseComptonProjToProjCase,
546  idxPostStep,n_order);
547  }
548  }
549  else {
550  if (theForcedInteractionForGamma) {
551  pmanager->AddProcess(theForcedInteractionForGamma);
552  n_order++;
553  pmanager->SetProcessOrdering(theForcedInteractionForGamma,
554  idxPostStep,n_order);
555  pmanager->SetProcessOrdering(theForcedInteractionForGamma,
556  idxAlongStep,n_order);
557  }
558  }
559  }
560 
561  if (particleName == "gamma") {
562  if (fUse_compton) {
563  pmanager->AddDiscreteProcess(theComptonScattering);
565  RegisterEmProcess(theComptonScattering,particle);
566  }
567  if (fUse_peeffect) {
568  pmanager->AddDiscreteProcess(thePEEffect);
570  RegisterEmProcess(thePEEffect,particle);
571  }
572  if (fUse_gamma_conversion) {
573  pmanager->AddDiscreteProcess(theGammaConversion);
574  }
575  }
576 
577  if (particleName == "e+" && fUse_gamma_conversion) {//positron
578  G4VProcess* theeplusMultipleScattering = new G4eMultipleScattering();
579  G4VProcess* theeplusIonisation = new G4eIonisation();
580  G4VProcess* theeplusBremsstrahlung = new G4eBremsstrahlung();
581  G4VProcess* theeplusAnnihilation = new G4eplusAnnihilation();
582 
583  // add processes
584  pmanager->AddProcess(theeplusMultipleScattering);
585  pmanager->AddProcess(theeplusIonisation);
586  pmanager->AddProcess(theeplusBremsstrahlung);
587  pmanager->AddProcess(theeplusAnnihilation);
588 
589  // set ordering for AtRestDoIt
590  pmanager->SetProcessOrderingToFirst(theeplusAnnihilation, idxAtRest);
591 
592  // set ordering for AlongStepDoIt
593  pmanager->SetProcessOrdering(theeplusMultipleScattering,
594  idxAlongStep,1);
595  pmanager->SetProcessOrdering(theeplusIonisation, idxAlongStep,2);
596  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxAlongStep,3);
597 
598  // set ordering for PostStepDoIt
599  pmanager->SetProcessOrdering(theeplusMultipleScattering,
600  idxPostStep,1);
601  pmanager->SetProcessOrdering(theeplusIonisation,idxPostStep,2);
602  pmanager->SetProcessOrdering(theeplusBremsstrahlung,idxPostStep,3);
603  pmanager->SetProcessOrdering(theeplusAnnihilation,idxPostStep,4);
604  }
605  if (particleName == "proton" && fUse_pionisation) {
606  if (fUse_ms && fUse_pionisation) pmanager->AddProcess(thepMS);
607 
608  if (fUse_pionisation){
609  pmanager->AddProcess(fPIonisation);
611  RegisterEnergyLossProcess(fPIonisation,particle);
612  }
613 
614  G4int n_order=0;
615  if (fUse_ms && fUse_pionisation) {
616  n_order++;
617  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
618  }
619 
620  if (fUse_pionisation) {
621  n_order++;
622  pmanager->SetProcessOrdering(fPIonisation,idxAlongStep,n_order);
623  }
624 
625  n_order=0;
626  if (fUse_ms && fUse_pionisation) {
627  n_order++;
628  pmanager->SetProcessOrdering(thepMS, idxPostStep,n_order);
629  }
630 
631  if (fUse_pionisation) {
632  n_order++;
633  pmanager->SetProcessOrdering(fPIonisation,idxPostStep,n_order);
634  }
635 
636  }
637 
638  if (particleName == "adj_proton" && fUse_pionisation) {
639  G4ContinuousGainOfEnergy* theContinuousGainOfEnergy =0;
640  if (fUse_pionisation ) {
641  theContinuousGainOfEnergy= new G4ContinuousGainOfEnergy();
642  theContinuousGainOfEnergy->SetLossFluctuations(
643  fUse_egain_fluctuation);
644  theContinuousGainOfEnergy->SetDirectEnergyLossProcess(fPIonisation);
645  theContinuousGainOfEnergy->SetDirectParticle(G4Proton::Proton());
646  pmanager->AddProcess(theContinuousGainOfEnergy);
647  }
648 
649  G4int n_order=0;
650  if (fUse_ms) {
651  n_order++;
652  pmanager->AddProcess(thepMS);
653  pmanager->SetProcessOrdering(thepMS, idxAlongStep,n_order);
654  }
655 
656  n_order++;
657  pmanager->SetProcessOrdering(theContinuousGainOfEnergy,idxAlongStep,
658  n_order);
659 
660  n_order++;
661  G4AdjointAlongStepWeightCorrection* theAlongStepWeightCorrection =
663  pmanager->AddProcess(theAlongStepWeightCorrection);
664  pmanager->SetProcessOrdering(theAlongStepWeightCorrection,
665  idxAlongStep,
666  n_order);
667  n_order=0;
668  if (fUse_pionisation) {
669  pmanager->AddProcess(thepInverseIonisationProjToProjCase);
670  n_order++;
671  pmanager->SetProcessOrdering(
672  thepInverseIonisationProjToProjCase,
673  idxPostStep,n_order);
674  }
675 
676  if (fUse_ms && fUse_pionisation) {
677  n_order++;
678  pmanager->SetProcessOrdering(thepMS,idxPostStep,n_order);
679  }
680  }
681  }
682 }
static G4AdjointSimManager * GetInstance()
static G4AdjointGamma * AdjointGamma()
void SetProcessManager(G4ProcessManager *aProcessManager)
void SetProcessOrderingToFirst(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void SetNewDisplacementFlag(G4bool)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4AdjointElectron * AdjointElectron()
void SetLowEnergyLimit(G4double aVal)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double aVal)
void RegisterAdjointParticle(G4ParticleDefinition *aPartDef)
void SetDirectEnergyLossProcess(G4VEnergyLossProcess *aProcess)
void SetLossFluctuations(G4bool val)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetDirectParticle(G4ParticleDefinition *p)
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 G4AdjointProton * AdjointProton()
G4ProcessManager * GetProcessManager() const
void ConsiderParticleAsPrimary(const G4String &particle_name)
void RegisterAdjointBremModel(G4VEmAdjointModel *aAdjointBremModel)
virtual void SetCSBiasingFactor(G4double aVal)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
G4eIonisation * fEminusIonisation
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
static G4AdjointCSManager * GetAdjointCSManager()
G4hIonisation * fPIonisation

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructGeneral ( )
protected

Definition at line 687 of file G4AdjointPhysicsList.cc.

688 {
689  // Add Decay Process
690  G4Decay* theDecayProcess = new G4Decay();
692  particleIterator->reset();
693  while( (*particleIterator)() ){
694  G4ParticleDefinition* particle = particleIterator->value();
695  G4ProcessManager* pmanager = particle->GetProcessManager();
696  if (theDecayProcess->IsApplicable(*particle)) {
697  pmanager ->AddProcess(theDecayProcess);
698  // set ordering for PostStepDoIt and AtRestDoIt
699  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
700  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
701  }
702  }
703 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4ProcessManager * GetProcessManager() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructLeptons ( )
protected

Definition at line 108 of file G4AdjointPhysicsList.cc.

109 {
110  // leptons
115 
120 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
static G4NeutrinoE * NeutrinoEDefinition()
Definition: G4NeutrinoE.cc:80
static G4AntiNeutrinoMu * AntiNeutrinoMuDefinition()
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
static G4AntiNeutrinoE * AntiNeutrinoEDefinition()
static G4NeutrinoMu * NeutrinoMuDefinition()
Definition: G4NeutrinoMu.cc:80

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructMesons ( )
protected

Definition at line 124 of file G4AdjointPhysicsList.cc.

125 {
126 // mesons
138 }
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4KaonZero * KaonZeroDefinition()
Definition: G4KaonZero.cc:99
static G4AntiKaonZero * AntiKaonZeroDefinition()
static G4KaonZeroShort * KaonZeroShortDefinition()
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:103
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4KaonZeroLong * KaonZeroLongDefinition()
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4EtaPrime * EtaPrimeDefinition()
Definition: G4EtaPrime.cc:106
static G4Eta * EtaDefinition()
Definition: G4Eta.cc:104

Here is the call graph for this function:

Here is the caller graph for this function:

void G4AdjointPhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 71 of file G4AdjointPhysicsList.cc.

72 {
73  // In this method, static member functions should be called
74  // for all particles which you want to use.
75  // This ensures that objects of these particle types will be
76  // created in the program.
82 }

Here is the call graph for this function:

void G4AdjointPhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 170 of file G4AdjointPhysicsList.cc.

Here is the call graph for this function:

void G4AdjointPhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 707 of file G4AdjointPhysicsList.cc.

708 {
709  if (verboseLevel >0){
710  G4cout << "G4AdjointPhysicsList::SetCuts:";
711  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
712  }
713 
714  // set cut values for gamma at first and for e- second and next for e+,
715  // because some processes for e+/e- need cut values for gamma
716  //
717  SetCutValue(defaultCutValue, "gamma");
720 
722 }
void SetCutValue(G4double aCut, const G4String &pname)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4AdjointPhysicsList::SetEmaxAdjModels ( G4double  aVal)
inline

Definition at line 75 of file G4AdjointPhysicsList.hh.

75 { fEmax_adj_models = aVal;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetEminAdjModels ( G4double  aVal)
inline

Definition at line 74 of file G4AdjointPhysicsList.hh.

74 { fEmin_adj_models = aVal;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetLossFluctuationFlag ( bool  aBool)

Definition at line 86 of file G4AdjointPhysicsList.cc.

87 {
89 }
void SetLossFluctuations(G4bool val)
G4eIonisation * fEminusIonisation

Here is the call graph for this function:

void G4AdjointPhysicsList::SetUseBrem ( bool  aBool)
inline

Definition at line 66 of file G4AdjointPhysicsList.hh.

66 {fUse_brem = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUseCompton ( bool  aBool)
inline

Definition at line 67 of file G4AdjointPhysicsList.hh.

67 {fUse_compton = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUseEgainFluctuation ( bool  aBool)
inline

Definition at line 72 of file G4AdjointPhysicsList.hh.

72  { fUse_egain_fluctuation
73  = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUseGammaConversion ( bool  aBool)
inline

Definition at line 70 of file G4AdjointPhysicsList.hh.

70  { fUse_gamma_conversion
71  = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUseIonisation ( bool  aBool)
inline

Definition at line 64 of file G4AdjointPhysicsList.hh.

64 {fUse_eionisation = aBool;}
void G4AdjointPhysicsList::SetUseMS ( bool  aBool)
inline

Definition at line 68 of file G4AdjointPhysicsList.hh.

68 {fUse_ms = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUsePEEffect ( bool  aBool)
inline

Definition at line 69 of file G4AdjointPhysicsList.hh.

69 {fUse_peeffect = aBool;}

Here is the caller graph for this function:

void G4AdjointPhysicsList::SetUseProtonIonisation ( bool  aBool)
inline

Definition at line 65 of file G4AdjointPhysicsList.hh.

65 {fUse_pionisation = aBool;}

Here is the caller graph for this function:

Member Data Documentation

G4eIonisation* G4AdjointPhysicsList::fEminusIonisation
protected

Definition at line 93 of file G4AdjointPhysicsList.hh.

G4hIonisation* G4AdjointPhysicsList::fPIonisation
protected

Definition at line 94 of file G4AdjointPhysicsList.hh.


The documentation for this class was generated from the following files: