Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
G4AdjointPhysicsList Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for G4AdjointPhysicsList:
G4VUserPhysicsList

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)
 

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 Attributes

G4eIonisationfEminusIonisation
 
G4hIonisationfPIonisation
 

Private Attributes

G4AdjointPhysicsMessengerfPhysicsMessenger
 
G4bool fUse_forced_interaction
 
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
 

Detailed Description

Definition at line 167 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ G4AdjointPhysicsList()

G4AdjointPhysicsList::G4AdjointPhysicsList ( )

Definition at line 49 of file G4AdjointPhysicsList.cc.

54 fUse_brem(true),fUse_compton(true),fUse_ms(true),
56 fEmin_adj_models(1.*keV), fEmax_adj_models(1.*MeV),
59{
60 defaultCutValue = 1.0*mm;
61 SetVerboseLevel(1);
63}
G4AdjointPhysicsMessenger * fPhysicsMessenger
G4eIonisation * fEminusIonisation

◆ ~G4AdjointPhysicsList()

G4AdjointPhysicsList::~G4AdjointPhysicsList ( )
virtual

Definition at line 67 of file G4AdjointPhysicsList.cc.

68{
69}

Member Function Documentation

◆ SetLossFluctuationFlag()

void G4AdjointPhysicsList::SetLossFluctuationFlag ( bool  aBool)

Definition at line 85 of file G4AdjointPhysicsList.cc.

86{
87 if (fEminusIonisation) fEminusIonisation->SetLossFluctuations(aBool);
88}

◆ SetUseIonisation()

void G4AdjointPhysicsList::SetUseIonisation ( bool  aBool)
inline

Definition at line 63 of file G4AdjointPhysicsList.hh.

63{fUse_eionisation = aBool;}

◆ SetUseProtonIonisation()

void G4AdjointPhysicsList::SetUseProtonIonisation ( bool  aBool)
inline

Definition at line 64 of file G4AdjointPhysicsList.hh.

64{fUse_pionisation = aBool;}

◆ SetUseBrem()

void G4AdjointPhysicsList::SetUseBrem ( bool  aBool)
inline

Definition at line 65 of file G4AdjointPhysicsList.hh.

65{fUse_brem = aBool;}

◆ SetUseCompton()

void G4AdjointPhysicsList::SetUseCompton ( bool  aBool)
inline

Definition at line 66 of file G4AdjointPhysicsList.hh.

66{fUse_compton = aBool;}

◆ SetUseMS()

void G4AdjointPhysicsList::SetUseMS ( bool  aBool)
inline

Definition at line 67 of file G4AdjointPhysicsList.hh.

67{fUse_ms = aBool;}

◆ SetUsePEEffect()

void G4AdjointPhysicsList::SetUsePEEffect ( bool  aBool)
inline

Definition at line 68 of file G4AdjointPhysicsList.hh.

68{fUse_peeffect = aBool;}

◆ SetUseGammaConversion()

void G4AdjointPhysicsList::SetUseGammaConversion ( bool  aBool)
inline

Definition at line 69 of file G4AdjointPhysicsList.hh.

◆ SetUseEgainFluctuation()

void G4AdjointPhysicsList::SetUseEgainFluctuation ( bool  aBool)
inline

Definition at line 71 of file G4AdjointPhysicsList.hh.

72 = aBool;}

◆ SetEminAdjModels()

void G4AdjointPhysicsList::SetEminAdjModels ( G4double  aVal)
inline

Definition at line 73 of file G4AdjointPhysicsList.hh.

73{ fEmin_adj_models = aVal;}

◆ SetEmaxAdjModels()

void G4AdjointPhysicsList::SetEmaxAdjModels ( G4double  aVal)
inline

Definition at line 74 of file G4AdjointPhysicsList.hh.

74{ fEmax_adj_models = aVal;}

◆ ConstructParticle()

void G4AdjointPhysicsList::ConstructParticle ( )
protectedvirtual

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}

◆ ConstructProcess()

void G4AdjointPhysicsList::ConstructProcess ( )
protectedvirtual

Definition at line 169 of file G4AdjointPhysicsList.cc.

170{
171 AddTransportation();
172 ConstructEM();
174}

◆ SetCuts()

void G4AdjointPhysicsList::SetCuts ( )
protectedvirtual

Definition at line 704 of file G4AdjointPhysicsList.cc.

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");
715 SetCutValue(defaultCutValue, "e-");
716 SetCutValue(defaultCutValue, "e+");
717
718 if (verboseLevel>0) DumpCutValuesTable();
719}

◆ ConstructBosons()

void G4AdjointPhysicsList::ConstructBosons ( )
protected

Definition at line 92 of file G4AdjointPhysicsList.cc.

93{
94 // pseudo-particles
95 G4Geantino::GeantinoDefinition();
96 G4ChargedGeantino::ChargedGeantinoDefinition();
97
98 // gamma
99 G4Gamma::GammaDefinition();
100
101 // optical photon
102 G4OpticalPhoton::OpticalPhotonDefinition();
103}

◆ ConstructLeptons()

void G4AdjointPhysicsList::ConstructLeptons ( )
protected

Definition at line 107 of file G4AdjointPhysicsList.cc.

108{
109 // leptons
110 G4Electron::ElectronDefinition();
111 G4Positron::PositronDefinition();
112 G4MuonPlus::MuonPlusDefinition();
113 G4MuonMinus::MuonMinusDefinition();
114
115 G4NeutrinoE::NeutrinoEDefinition();
116 G4AntiNeutrinoE::AntiNeutrinoEDefinition();
117 G4NeutrinoMu::NeutrinoMuDefinition();
118 G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
119}

◆ ConstructMesons()

void G4AdjointPhysicsList::ConstructMesons ( )
protected

Definition at line 123 of file G4AdjointPhysicsList.cc.

124{
125// mesons
126 G4PionPlus::PionPlusDefinition();
127 G4PionMinus::PionMinusDefinition();
128 G4PionZero::PionZeroDefinition();
129 G4Eta::EtaDefinition();
130 G4EtaPrime::EtaPrimeDefinition();
131 G4KaonPlus::KaonPlusDefinition();
132 G4KaonMinus::KaonMinusDefinition();
133 G4KaonZero::KaonZeroDefinition();
134 G4AntiKaonZero::AntiKaonZeroDefinition();
135 G4KaonZeroLong::KaonZeroLongDefinition();
136 G4KaonZeroShort::KaonZeroShortDefinition();
137}

◆ ConstructBaryons()

void G4AdjointPhysicsList::ConstructBaryons ( )
protected

Definition at line 141 of file G4AdjointPhysicsList.cc.

142{
143// barions
144 G4Proton::ProtonDefinition();
145 G4AntiProton::AntiProtonDefinition();
146 G4Neutron::NeutronDefinition();
147 G4AntiNeutron::AntiNeutronDefinition();
148}

◆ ConstructAdjointParticles()

void G4AdjointPhysicsList::ConstructAdjointParticles ( )
protected

Definition at line 155 of file G4AdjointPhysicsList.cc.

156{
157// adjoint_gammma
158 G4AdjointGamma::AdjointGammaDefinition();
159
160// adjoint_electron
161 G4AdjointElectron::AdjointElectronDefinition();
162
163// adjoint_proton
164 G4AdjointProton::AdjointProtonDefinition();
165}

◆ ConstructGeneral()

void G4AdjointPhysicsList::ConstructGeneral ( )
protected

Definition at line 684 of file G4AdjointPhysicsList.cc.

685{
686 // Add Decay Process
687 G4Decay* theDecayProcess = new G4Decay();
688 auto particleIterator=GetParticleIterator();
689 particleIterator->reset();
690 while( (*particleIterator)() ){
691 G4ParticleDefinition* particle = particleIterator->value();
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}

◆ ConstructEM()

void G4AdjointPhysicsList::ConstructEM ( )
protected

Definition at line 222 of file G4AdjointPhysicsList.cc.

223{
224 G4AdjointCSManager* theCSManager =
225 G4AdjointCSManager::GetAdjointCSManager();
226 G4AdjointSimManager* theAdjointSimManager =
227 G4AdjointSimManager::GetInstance();
228
229 theCSManager->RegisterAdjointParticle(
230 G4AdjointElectron::AdjointElectron());
231
233 theCSManager->RegisterAdjointParticle(
234 G4AdjointGamma::AdjointGamma());
235
236 if (fUse_eionisation) {
238 new G4eIonisation();
239 fEminusIonisation->SetLossFluctuations(fUse_egain_fluctuation);
240 }
241 if (fUse_pionisation) {
243 fPIonisation->SetLossFluctuations(fUse_egain_fluctuation);
244 theCSManager->RegisterAdjointParticle(
245 G4AdjointProton::AdjointProton());
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;
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 //--------------------------------------------------------
396 auto particleIterator=GetParticleIterator();
397 particleIterator->reset();
398 while( (*particleIterator)() ){
399 G4ParticleDefinition* particle = particleIterator->value();
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);
411 G4AdjointCSManager::GetAdjointCSManager()->
412 RegisterEnergyLossProcess(fEminusIonisation,particle);
413 }
415 pmanager->AddProcess(theeminusBremsstrahlung);
416 G4AdjointCSManager::GetAdjointCSManager()->
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++;
426 pmanager->SetProcessOrdering(fEminusIonisation,idxAlongStep,n_order);
427 }
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++;
440 pmanager->SetProcessOrdering(fEminusIonisation,idxPostStep,n_order);
441 }
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 }
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
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);
561 G4AdjointCSManager::GetAdjointCSManager()->
562 RegisterEmProcess(theComptonScattering,particle);
563 }
564 if (fUse_peeffect) {
565 pmanager->AddDiscreteProcess(thePEEffect);
566 G4AdjointCSManager::GetAdjointCSManager()->
567 RegisterEmProcess(thePEEffect,particle);
568 }
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);
607 G4AdjointCSManager::GetAdjointCSManager()->
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}

Member Data Documentation

◆ fEminusIonisation

G4eIonisation* G4AdjointPhysicsList::fEminusIonisation
protected

Definition at line 92 of file G4AdjointPhysicsList.hh.

◆ fPIonisation

G4hIonisation* G4AdjointPhysicsList::fPIonisation
protected

Definition at line 93 of file G4AdjointPhysicsList.hh.

◆ fPhysicsMessenger

G4AdjointPhysicsMessenger* G4AdjointPhysicsList::fPhysicsMessenger
private

Definition at line 96 of file G4AdjointPhysicsList.hh.

◆ fUse_forced_interaction

G4bool G4AdjointPhysicsList::fUse_forced_interaction
private

Definition at line 97 of file G4AdjointPhysicsList.hh.

◆ fUse_eionisation

G4bool G4AdjointPhysicsList::fUse_eionisation
private

Definition at line 98 of file G4AdjointPhysicsList.hh.

◆ fUse_pionisation

G4bool G4AdjointPhysicsList::fUse_pionisation
private

Definition at line 99 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_ms

G4bool G4AdjointPhysicsList::fUse_ms
private

Definition at line 102 of file G4AdjointPhysicsList.hh.

◆ fUse_egain_fluctuation

G4bool G4AdjointPhysicsList::fUse_egain_fluctuation
private

Definition at line 103 of file G4AdjointPhysicsList.hh.

◆ fUse_peeffect

G4bool G4AdjointPhysicsList::fUse_peeffect
private

Definition at line 104 of file G4AdjointPhysicsList.hh.

◆ fUse_gamma_conversion

G4bool G4AdjointPhysicsList::fUse_gamma_conversion
private

Definition at line 105 of file G4AdjointPhysicsList.hh.

◆ fEmin_adj_models

G4double G4AdjointPhysicsList::fEmin_adj_models
private

Definition at line 106 of file G4AdjointPhysicsList.hh.

◆ fEmax_adj_models

G4double G4AdjointPhysicsList::fEmax_adj_models
private

Definition at line 107 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_brem

G4double G4AdjointPhysicsList::fCS_biasing_factor_brem
private

Definition at line 109 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.


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

Applications | User Support | Publications | Collaboration