Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
GarfieldMessenger Class Reference

#include <GarfieldMessenger.hh>

+ Inheritance diagram for GarfieldMessenger:

Public Member Functions

 GarfieldMessenger (GarfieldDetectorConstruction *)
 
 ~GarfieldMessenger ()
 
virtual void SetNewValue (G4UIcommand *, G4String)
 

Detailed Description

Definition at line 46 of file GarfieldMessenger.hh.

Constructor & Destructor Documentation

◆ GarfieldMessenger()

GarfieldMessenger::GarfieldMessenger ( GarfieldDetectorConstruction * Det)

Definition at line 40 of file GarfieldMessenger.cc.

41 : G4UImessenger(),
42 fDetector(Det),
43 fExampleDir(0),
44 fAbsorberDir(0),
45 fGarfieldPhysicsDir(0),
46 fMaterialCmd(0),
47 fIsotopeCmd(0),
48 fIonizationModelCmd(0),
49 fGeant4ParticleCmd(0) {
50 fExampleDir = new G4UIdirectory("/exampleGarfield/");
51 fExampleDir->SetGuidance("Commands specific to this example");
52
53 G4bool broadcast = false;
54 fAbsorberDir = new G4UIdirectory("/exampleGarfield/absorber/", broadcast);
55 fAbsorberDir->SetGuidance("Absorber construction commands");
56
57 fGarfieldPhysicsDir =
58 new G4UIdirectory("/exampleGarfield/physics/", broadcast);
59 fGarfieldPhysicsDir->SetGuidance(
60 "Particle and energy ranges for Garfield++ physics model");
61
62 fMaterialCmd =
63 new G4UIcmdWithAString("/exampleGarfield/absorber/setMat", this);
64 fMaterialCmd->SetGuidance("Select material of the absorber:");
65 fMaterialCmd->SetParameterName("choice", false);
66 fMaterialCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
67
68 fIsotopeCmd =
69 new G4UIcommand("/exampleGarfield/absorber/setIsotopeMat", this);
70 fIsotopeCmd->SetGuidance("Build and select a material with single isotope");
71 fIsotopeCmd->SetGuidance(" symbol of isotope, Z, A, density of material");
72 //
73 G4UIparameter* symbPrm = new G4UIparameter("isotope", 's', false);
74 symbPrm->SetGuidance("isotope symbol");
75 fIsotopeCmd->SetParameter(symbPrm);
76 //
77 G4UIparameter* ZPrm = new G4UIparameter("Z", 'i', false);
78 ZPrm->SetGuidance("Z");
79 ZPrm->SetParameterRange("Z>0");
80 fIsotopeCmd->SetParameter(ZPrm);
81 //
82 G4UIparameter* APrm = new G4UIparameter("A", 'i', false);
83 APrm->SetGuidance("A");
84 APrm->SetParameterRange("A>0");
85 fIsotopeCmd->SetParameter(APrm);
86 //
87 G4UIparameter* densityPrm = new G4UIparameter("density", 'd', false);
88 densityPrm->SetGuidance("density of material");
89 densityPrm->SetParameterRange("density>0.");
90 fIsotopeCmd->SetParameter(densityPrm);
91 //
92 G4UIparameter* unitPrm = new G4UIparameter("unit", 's', false);
93 unitPrm->SetGuidance("unit of density");
94 G4String unitList = G4UIcommand::UnitsList(G4UIcommand::CategoryOf("g/cm3"));
95 unitPrm->SetParameterCandidates(unitList);
96 fIsotopeCmd->SetParameter(unitPrm);
97 //
98 fIsotopeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
99
100 fIonizationModelCmd =
101 new G4UIcommand("/exampleGarfield/physics/setIonizationModel", this);
102 fIonizationModelCmd->SetGuidance("Select ionization model for Garfield++");
103 fIonizationModelCmd->SetGuidance(
104 " and choose whether to use default particles");
105 fIonizationModelCmd->SetGuidance(" and energy ranges for the chosen model");
106 //
107 G4UIparameter* ionizationModelPrm =
108 new G4UIparameter("ionizationModel", 's', false);
109 ionizationModelPrm->SetGuidance(
110 "ionization model (1. PAIPhot, 2. PAI or 3. Heed)");
111 ionizationModelPrm->SetGuidance(
112 " 1. PAIPhot model in Geant4, delta electrons transported by Heed");
113 ionizationModelPrm->SetGuidance(
114 " 2. PAI model in Geant4, delta electrons transported by Heed");
115 ionizationModelPrm->SetGuidance(" 3. Use directly Heed");
116 fIonizationModelCmd->SetParameter(ionizationModelPrm);
117 //
118 G4UIparameter* useDefaultsPrm = new G4UIparameter("useDefaults", 'b', false);
119 useDefaultsPrm->SetGuidance(
120 "true to use default, false to manually choose particles and energies");
121 fIonizationModelCmd->SetParameter(useDefaultsPrm);
122 //
123 fIonizationModelCmd->AvailableForStates(G4State_PreInit);
124
125 fGeant4ParticleCmd = new G4UIcommand(
126 "/exampleGarfield/physics/setGeant4ParticleTypeAndEnergy", this);
127 fGeant4ParticleCmd->SetGuidance(
128 "Select particle types and energies for PAI and PAIPhot model");
129 fGeant4ParticleCmd->SetGuidance(" in Geant4");
130 //
131 G4UIparameter* particleGeant4Prm =
132 new G4UIparameter("particleName", 's', false);
133 particleGeant4Prm->SetGuidance(
134 "Particle name (e-, e+, mu-, mu+, proton, pi-, pi+, alpha, He3, "
135 "GenericIon)");
136 fGeant4ParticleCmd->SetParameter(particleGeant4Prm);
137 //
138 G4UIparameter* minEnergyGeant4Prm =
139 new G4UIparameter("minimumEnergyGeant4", 'd', false);
140 minEnergyGeant4Prm->SetGuidance("minimum energy");
141 minEnergyGeant4Prm->SetParameterRange("minimumEnergyGeant4>=0");
142 fGeant4ParticleCmd->SetParameter(minEnergyGeant4Prm);
143 //
144 G4UIparameter* maxEnergyGeant4Prm =
145 new G4UIparameter("maximumEnergyGeant4", 'd', false);
146 maxEnergyGeant4Prm->SetGuidance("maximum energy");
147 maxEnergyGeant4Prm->SetParameterRange("maximumEnergyGeant4>=0");
148 fGeant4ParticleCmd->SetParameter(maxEnergyGeant4Prm);
149 //
150 G4UIparameter* unitGeant4Prm = new G4UIparameter("unit", 's', false);
151 unitGeant4Prm->SetGuidance("unit of energy");
152 G4String unitListGeant4 =
153 G4UIcommand::UnitsList(G4UIcommand::CategoryOf("MeV"));
154 unitGeant4Prm->SetParameterCandidates(unitListGeant4);
155 fGeant4ParticleCmd->SetParameter(unitGeant4Prm);
156 //
157 fGeant4ParticleCmd->AvailableForStates(G4State_PreInit);
158
159 fGarfieldParticleCmd = new G4UIcommand(
160 "/exampleGarfield/physics/setGarfieldParticleTypeAndEnergy", this);
161 fGarfieldParticleCmd->SetGuidance(
162 "Select particle types and energies for Heed model.");
163 fGarfieldParticleCmd->SetGuidance(
164 " For PAI and PAIPhot model choose at which energy electrons are");
165 fGarfieldParticleCmd->SetGuidance(
166 " transported as delta electrons by Heed, and treatment of gammas");
167 //
168 G4UIparameter* particleGarfieldPrm =
169 new G4UIparameter("particleName", 's', false);
170 particleGarfieldPrm->SetGuidance(
171 "Particle name (gamma, e-, e+, mu-, mu+, proton, anti_proton, pi-, pi+, "
172 "kaon, kaon+, alpha, deuteron)");
173 fGarfieldParticleCmd->SetParameter(particleGarfieldPrm);
174 //
175 G4UIparameter* minEnergyGarfieldPrm =
176 new G4UIparameter("minimumEnergyGarfield", 'd', false);
177 minEnergyGarfieldPrm->SetGuidance("minimum energy");
178 minEnergyGarfieldPrm->SetParameterRange("minimumEnergyGarfield>=0");
179 fGarfieldParticleCmd->SetParameter(minEnergyGarfieldPrm);
180 //
181 G4UIparameter* maxEnergyGarfieldPrm =
182 new G4UIparameter("maximumEnergyGarfield", 'd', false);
183 maxEnergyGarfieldPrm->SetGuidance("maximum energy");
184 maxEnergyGarfieldPrm->SetParameterRange("maximumEnergyGarfield>=0");
185 fGarfieldParticleCmd->SetParameter(maxEnergyGarfieldPrm);
186 //
187 G4UIparameter* unitGarfieldPrm = new G4UIparameter("unit", 's', false);
188 unitGarfieldPrm->SetGuidance("unit of energy");
189 G4String unitListGarfield =
190 G4UIcommand::UnitsList(G4UIcommand::CategoryOf("MeV"));
191 unitGarfieldPrm->SetParameterCandidates(unitListGarfield);
192 fGarfieldParticleCmd->SetParameter(unitGarfieldPrm);
193 //
194 fGarfieldParticleCmd->AvailableForStates(G4State_PreInit);
195}

◆ ~GarfieldMessenger()

GarfieldMessenger::~GarfieldMessenger ( )

Definition at line 199 of file GarfieldMessenger.cc.

199 {
200 delete fMaterialCmd;
201 delete fIsotopeCmd;
202 delete fAbsorberDir;
203 delete fExampleDir;
204 delete fGarfieldPhysicsDir;
205 delete fIonizationModelCmd;
206 delete fGeant4ParticleCmd;
207 delete fGarfieldParticleCmd;
208}

Member Function Documentation

◆ SetNewValue()

void GarfieldMessenger::SetNewValue ( G4UIcommand * command,
G4String newValue )
virtual

Definition at line 212 of file GarfieldMessenger.cc.

212 {
213 if (command == fMaterialCmd) {
214 fDetector->SetAbsorberMaterial(newValue);
215 } else if (command == fIsotopeCmd) {
216 G4int Z;
217 G4int A;
218 G4double dens;
219 G4String name, unt;
220 std::istringstream is(newValue);
221 is >> name >> Z >> A >> dens >> unt;
222 dens *= G4UIcommand::ValueOf(unt);
223 fDetector->AbsorberMaterialWithSingleIsotope(name, name, dens, Z, A);
224 fDetector->SetAbsorberMaterial(name);
225 } else if (command == fIonizationModelCmd) {
226 GarfieldPhysics* garfieldPhysics = GarfieldPhysics::GetInstance();
227 G4String modelName;
228 G4bool useDefaults;
229 std::istringstream is(newValue);
230 is >> modelName >> std::boolalpha >> useDefaults;
231 garfieldPhysics->SetIonizationModel(modelName, useDefaults);
232 } else if (command == fGeant4ParticleCmd) {
233 GarfieldPhysics* garfieldPhysics = GarfieldPhysics::GetInstance();
234 G4String particleName, unit, programName;
235 G4double minEnergy;
236 G4double maxEnergy;
237 std::istringstream is(newValue);
238 is >> particleName >> minEnergy >> maxEnergy >> unit;
239 minEnergy *= G4UIcommand::ValueOf(unit);
240 maxEnergy *= G4UIcommand::ValueOf(unit);
241 garfieldPhysics->AddParticleName(particleName, minEnergy, maxEnergy,
242 "geant4");
243 } else if (command == fGarfieldParticleCmd) {
244 GarfieldPhysics* garfieldPhysics = GarfieldPhysics::GetInstance();
245 G4String particleName, unit, programName;
246 G4double minEnergy;
247 G4double maxEnergy;
248 std::istringstream is(newValue);
249 is >> particleName >> minEnergy >> maxEnergy >> unit;
250 minEnergy *= G4UIcommand::ValueOf(unit);
251 maxEnergy *= G4UIcommand::ValueOf(unit);
252 garfieldPhysics->AddParticleName(particleName, minEnergy, maxEnergy,
253 "garfield");
254 }
255}
static GarfieldPhysics * GetInstance()
void SetIonizationModel(std::string model, bool useDefaults=true)
void AddParticleName(const std::string particleName, double ekin_min_MeV, double ekin_max_MeV, std::string program)

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