Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutrinoPhysicsMessenger Class Reference

#include <G4NeutrinoPhysicsMessenger.hh>

Inheritance diagram for G4NeutrinoPhysicsMessenger:

Public Member Functions

 G4NeutrinoPhysicsMessenger (G4NeutrinoPhysics *af)
 ~G4NeutrinoPhysicsMessenger () override
void SetNewValue (G4UIcommand *aComm, G4String aS) override
G4NeutrinoPhysicsMessengeroperator= (const G4NeutrinoPhysicsMessenger &right)=delete
 G4NeutrinoPhysicsMessenger (const G4NeutrinoPhysicsMessenger &)=delete
Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
virtual ~G4UImessenger ()
virtual G4String GetCurrentValue (G4UIcommand *command)
G4bool CommandsShouldBeInMaster () const

Additional Inherited Members

Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
G4String LtoS (G4long l)
G4String DtoS (G4double a)
G4String BtoS (G4bool b)
G4int StoI (const G4String &s)
G4long StoL (const G4String &s)
G4double StoD (const G4String &s)
G4bool StoB (const G4String &s)
void AddUIcommand (G4UIcommand *newCommand)
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
template<typename T>
T * CreateCommand (const G4String &cname, const G4String &dsc)
Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
G4String baseDirName = ""
G4bool commandsShouldBeInMaster = false

Detailed Description

Definition at line 49 of file G4NeutrinoPhysicsMessenger.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoPhysicsMessenger() [1/2]

G4NeutrinoPhysicsMessenger::G4NeutrinoPhysicsMessenger ( G4NeutrinoPhysics * af)
explicit

Definition at line 41 of file G4NeutrinoPhysicsMessenger.cc.

42 : theB(ab)
43{
44 // general stuff.
45 aDir = new G4UIdirectory("/physics_lists/nu/", false);
46 aDir->SetGuidance("tailoring the neutrino processes.");
47
48 theNu = new G4UIcmdWithABool("/physics_lists/nu/NeutrinoActivation",this);
49 theNu->SetGuidance("Activation of neutrino-nucleus processes");
50 theNu->AvailableForStates(G4State_PreInit);
51 theNu->SetToBeBroadcasted(false);
52
53 theNuETX = new G4UIcmdWithABool("/physics_lists/nu/NuETotXscActivation",this);
54 theNuETX->SetGuidance("Activation of neutrino-electron processes");
55 theNuETX->AvailableForStates(G4State_PreInit);
56 theNuETX->SetToBeBroadcasted(false);
57
58 theNuEleCcBF = new G4UIcmdWithADouble("/physics_lists/nu/NuEleCcBias",this);
59 theNuEleCcBF->SetGuidance("Neutrino-electron charge current bias factor");
60 theNuEleCcBF->AvailableForStates(G4State_PreInit);
61 theNuEleCcBF->SetToBeBroadcasted(false);
62
63 theNuEleNcBF = new G4UIcmdWithADouble("/physics_lists/nu/NuEleNcBias",this);
64 theNuEleNcBF->SetGuidance("Neutrino-electron neutral current bias factor");
65 theNuEleNcBF->AvailableForStates(G4State_PreInit);
66 theNuEleNcBF->SetToBeBroadcasted(false);
67
68 theNuNucleusBF = new G4UIcmdWithADouble("/physics_lists/nu/NuNucleusBias",this);
69 theNuNucleusBF->SetGuidance("Neutrino-nucleus cross section bias factor");
70 theNuNucleusBF->AvailableForStates(G4State_PreInit);
71 theNuNucleusBF->SetToBeBroadcasted(false);
72
73 theNuOscDistanceBF = new G4UIcmdWithADouble("/physics_lists/nu/NuOscDistanceBias",this);
74 theNuOscDistanceBF->SetGuidance("Neutrino-oscillation distance bias factor");
75 theNuOscDistanceBF->AvailableForStates(G4State_PreInit);
76 theNuOscDistanceBF->SetToBeBroadcasted(false);
77
78 theNuDN = new G4UIcmdWithAString("/physics_lists/nu/NuDetectorName",this);
79 theNuDN->SetGuidance("Set neutrino detector name");
80 theNuDN->AvailableForStates(G4State_PreInit);
81 theNuDN->SetToBeBroadcasted(false);
82
83 theNuODN = new G4UIcmdWithAString("/physics_lists/nu/NuOscDistanceName",this);
84 theNuODN->SetGuidance("Set neutrino oscillation distance region name");
85 theNuODN->AvailableForStates(G4State_PreInit);
86 theNuODN->SetToBeBroadcasted(false);
87}
@ G4State_PreInit

Referenced by G4NeutrinoPhysicsMessenger(), and operator=().

◆ ~G4NeutrinoPhysicsMessenger()

G4NeutrinoPhysicsMessenger::~G4NeutrinoPhysicsMessenger ( )
override

Definition at line 89 of file G4NeutrinoPhysicsMessenger.cc.

90{
91 delete theNu;
92 delete theNuETX;
93
94 delete theNuEleCcBF;
95 delete theNuEleNcBF;
96 delete theNuNucleusBF;
97 delete theNuOscDistanceBF;
98
99 delete theNuDN;
100 delete theNuODN;
101
102 delete aDir;
103}

◆ G4NeutrinoPhysicsMessenger() [2/2]

G4NeutrinoPhysicsMessenger::G4NeutrinoPhysicsMessenger ( const G4NeutrinoPhysicsMessenger & )
delete

Member Function Documentation

◆ operator=()

G4NeutrinoPhysicsMessenger & G4NeutrinoPhysicsMessenger::operator= ( const G4NeutrinoPhysicsMessenger & right)
delete

◆ SetNewValue()

void G4NeutrinoPhysicsMessenger::SetNewValue ( G4UIcommand * aComm,
G4String aS )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 105 of file G4NeutrinoPhysicsMessenger.cc.

106{
107 if (aComm==theNuETX)
108 theB->NuETotXscActivated(theNuETX->GetNewBoolValue(aS));
109 else if (aComm==theNuEleCcBF)
110 theB->SetNuEleCcBias(theNuEleCcBF->GetNewDoubleValue(aS));
111 else if (aComm==theNuEleNcBF)
112 theB->SetNuEleNcBias(theNuEleNcBF->GetNewDoubleValue(aS));
113 else if (aComm==theNuNucleusBF)
114 theB->SetNuNucleusBias(theNuNucleusBF->GetNewDoubleValue(aS));
115 else if (aComm==theNuOscDistanceBF)
116 theB->SetNuOscDistanceBias(theNuOscDistanceBF->GetNewDoubleValue(aS));
117 else if(aComm==theNuDN)
118 theB->SetNuDetectorName(aS);
119 else if(aComm==theNuODN)
120 theB->SetNuOscDistanceName(aS);
121}

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