Geant4  10.02.p03
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
 

Private Attributes

G4AdjointPhysicsMessengerfPhysicsMessenger
 
G4bool fUse_eionisation
 
G4bool fUse_pionisation
 
G4bool fUse_brem
 
G4bool fUse_compton
 
G4bool fUse_ms
 
G4bool fUse_egain_fluctuation
 
G4bool fUse_peeffect
 
G4bool fUse_gamma_conversion
 
G4double fEmin_adj_models
 
G4double fEmax_adj_models
 
G4double fCS_biasing_factor_compton
 
G4double fCS_biasing_factor_brem
 
G4double fCS_biasing_factor_ionisation
 
G4double fCS_biasing_factor_PEeffect
 

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::G4AdjointPhysicsList ( )

Definition at line 50 of file G4AdjointPhysicsList.cc.

54  fUse_brem(true),fUse_compton(true),fUse_ms(true),
59 {
60  defaultCutValue = 1.0*mm;
61  SetVerboseLevel(1);
63 }
static const double MeV
Definition: G4SIunits.hh:211
G4AdjointPhysicsMessenger * fPhysicsMessenger
void SetVerboseLevel(G4int value)
G4eIonisation * fEminusIonisation
static const double keV
Definition: G4SIunits.hh:213
static const double mm
Definition: G4SIunits.hh:114
G4hIonisation * fPIonisation
Here is the call graph for this function:

◆ ~G4AdjointPhysicsList()

G4AdjointPhysicsList::~G4AdjointPhysicsList ( )
virtual

Definition at line 67 of file G4AdjointPhysicsList.cc.

68 {
69 }

Member Function Documentation

◆ ConstructAdjointParticles()

void G4AdjointPhysicsList::ConstructAdjointParticles ( )
protected

Definition at line 155 of file G4AdjointPhysicsList.cc.

156 {
157 // adjoint_gammma
159 
160 // adjoint_electron
162 
163 // adjoint_proton
165 }
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:

◆ ConstructBaryons()

void G4AdjointPhysicsList::ConstructBaryons ( )
protected

Definition at line 141 of file G4AdjointPhysicsList.cc.

142 {
143 // barions
148 }
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:

◆ ConstructBosons()

void G4AdjointPhysicsList::ConstructBosons ( )
protected

Definition at line 92 of file G4AdjointPhysicsList.cc.

93 {
94  // pseudo-particles
97 
98  // gamma
100 
101  // optical photon
103 }
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:

◆ ConstructEM()

void G4AdjointPhysicsList::ConstructEM ( )
protected

Definition at line 211 of file G4AdjointPhysicsList.cc.

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

◆ ConstructGeneral()

void G4AdjointPhysicsList::ConstructGeneral ( )
protected

Definition at line 641 of file G4AdjointPhysicsList.cc.

642 {
643  // Add Decay Process
644  G4Decay* theDecayProcess = new G4Decay();
646  particleIterator->reset();
647  while( (*particleIterator)() ){
648  G4ParticleDefinition* particle = particleIterator->value();
649  G4ProcessManager* pmanager = particle->GetProcessManager();
650  if (theDecayProcess->IsApplicable(*particle)) {
651  pmanager ->AddProcess(theDecayProcess);
652  // set ordering for PostStepDoIt and AtRestDoIt
653  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
654  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
655  }
656  }
657 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4ProcessManager * GetProcessManager() const
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructLeptons()

void G4AdjointPhysicsList::ConstructLeptons ( )
protected

Definition at line 107 of file G4AdjointPhysicsList.cc.

108 {
109  // leptons
114 
119 }
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:

◆ ConstructMesons()

void G4AdjointPhysicsList::ConstructMesons ( )
protected

Definition at line 123 of file G4AdjointPhysicsList.cc.

124 {
125 // mesons
137 }
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:

◆ ConstructParticle()

void G4AdjointPhysicsList::ConstructParticle ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 70 of file G4AdjointPhysicsList.cc.

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 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructProcess()

void G4AdjointPhysicsList::ConstructProcess ( void  )
protectedvirtual

Implements G4VUserPhysicsList.

Definition at line 169 of file G4AdjointPhysicsList.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetCuts()

void G4AdjointPhysicsList::SetCuts ( )
protectedvirtual

Reimplemented from G4VUserPhysicsList.

Definition at line 661 of file G4AdjointPhysicsList.cc.

662 {
663  if (verboseLevel >0){
664  G4cout << "G4AdjointPhysicsList::SetCuts:";
665  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
666  }
667 
668  // set cut values for gamma at first and for e- second and next for e+,
669  // because some processes for e+/e- need cut values for gamma
670  //
671  SetCutValue(defaultCutValue, "gamma");
674 
676 }
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:
Here is the caller graph for this function:

◆ SetEmaxAdjModels()

void G4AdjointPhysicsList::SetEmaxAdjModels ( G4double  aVal)
inline

Definition at line 75 of file G4AdjointPhysicsList.hh.

75 { fEmax_adj_models = aVal;}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetEminAdjModels()

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:

◆ SetLossFluctuationFlag()

void G4AdjointPhysicsList::SetLossFluctuationFlag ( bool  aBool)

Definition at line 85 of file G4AdjointPhysicsList.cc.

86 {
88 }
void SetLossFluctuations(G4bool val)
G4eIonisation * fEminusIonisation
Here is the call graph for this function:

◆ SetUseBrem()

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:

◆ SetUseCompton()

void G4AdjointPhysicsList::SetUseCompton ( bool  aBool)
inline

Definition at line 67 of file G4AdjointPhysicsList.hh.

Here is the caller graph for this function:

◆ SetUseEgainFluctuation()

void G4AdjointPhysicsList::SetUseEgainFluctuation ( bool  aBool)
inline

Definition at line 72 of file G4AdjointPhysicsList.hh.

73  = aBool;}
Here is the caller graph for this function:

◆ SetUseGammaConversion()

void G4AdjointPhysicsList::SetUseGammaConversion ( bool  aBool)
inline

Definition at line 70 of file G4AdjointPhysicsList.hh.

71  = aBool;}
Here is the caller graph for this function:

◆ SetUseIonisation()

void G4AdjointPhysicsList::SetUseIonisation ( bool  aBool)
inline

Definition at line 64 of file G4AdjointPhysicsList.hh.

◆ SetUseMS()

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:

◆ SetUsePEEffect()

void G4AdjointPhysicsList::SetUsePEEffect ( bool  aBool)
inline

Definition at line 69 of file G4AdjointPhysicsList.hh.

Here is the caller graph for this function:

◆ SetUseProtonIonisation()

void G4AdjointPhysicsList::SetUseProtonIonisation ( bool  aBool)
inline

Definition at line 65 of file G4AdjointPhysicsList.hh.

Here is the caller graph for this function:

Member Data Documentation

◆ fCS_biasing_factor_brem

G4double G4AdjointPhysicsList::fCS_biasing_factor_brem
private

Definition at line 109 of file G4AdjointPhysicsList.hh.

◆ fCS_biasing_factor_compton

G4double G4AdjointPhysicsList::fCS_biasing_factor_compton
private

Definition at line 108 of file G4AdjointPhysicsList.hh.

◆ fCS_biasing_factor_ionisation

G4double G4AdjointPhysicsList::fCS_biasing_factor_ionisation
private

Definition at line 110 of file G4AdjointPhysicsList.hh.

◆ fCS_biasing_factor_PEeffect

G4double G4AdjointPhysicsList::fCS_biasing_factor_PEeffect
private

Definition at line 111 of file G4AdjointPhysicsList.hh.

◆ fEmax_adj_models

G4double G4AdjointPhysicsList::fEmax_adj_models
private

Definition at line 107 of file G4AdjointPhysicsList.hh.

◆ fEmin_adj_models

G4double G4AdjointPhysicsList::fEmin_adj_models
private

Definition at line 106 of file G4AdjointPhysicsList.hh.

◆ fEminusIonisation

G4eIonisation* G4AdjointPhysicsList::fEminusIonisation
protected

Definition at line 93 of file G4AdjointPhysicsList.hh.

◆ fPhysicsMessenger

G4AdjointPhysicsMessenger* G4AdjointPhysicsList::fPhysicsMessenger
private

Definition at line 97 of file G4AdjointPhysicsList.hh.

◆ fPIonisation

G4hIonisation* G4AdjointPhysicsList::fPIonisation
protected

Definition at line 94 of file G4AdjointPhysicsList.hh.

◆ fUse_brem

G4bool G4AdjointPhysicsList::fUse_brem
private

Definition at line 100 of file G4AdjointPhysicsList.hh.

◆ fUse_compton

G4bool G4AdjointPhysicsList::fUse_compton
private

Definition at line 101 of file G4AdjointPhysicsList.hh.

◆ fUse_egain_fluctuation

G4bool G4AdjointPhysicsList::fUse_egain_fluctuation
private

Definition at line 103 of file G4AdjointPhysicsList.hh.

◆ fUse_eionisation

G4bool G4AdjointPhysicsList::fUse_eionisation
private

Definition at line 98 of file G4AdjointPhysicsList.hh.

◆ fUse_gamma_conversion

G4bool G4AdjointPhysicsList::fUse_gamma_conversion
private

Definition at line 105 of file G4AdjointPhysicsList.hh.

◆ fUse_ms

G4bool G4AdjointPhysicsList::fUse_ms
private

Definition at line 102 of file G4AdjointPhysicsList.hh.

◆ fUse_peeffect

G4bool G4AdjointPhysicsList::fUse_peeffect
private

Definition at line 104 of file G4AdjointPhysicsList.hh.

◆ fUse_pionisation

G4bool G4AdjointPhysicsList::fUse_pionisation
private

Definition at line 99 of file G4AdjointPhysicsList.hh.


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