VMC Examples Version 6.6
Loading...
Searching...
No Matches
GarfieldMessenger.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Virtual Monte Carlo examples
3// Copyright (C) 2007 - 2016 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
10/// \file ExGarfield/geant4/src/GarfieldG4FastSimulationModel.cxx
11/// \brief Implementation of the GarfieldG4FastSimulationModel class
12///
13/// Garfield garfieldpp example adapted to Virtual Monte Carlo.
14/// This class is imported from garfieldpp example.
15/// There is disabled WriteGeometryToGDML function in order to avoid
16/// requiring Geant4 optional library.
17///
18/// \date 28/10/2015
19/// \author D. Pheiffer, CERN
20
21#include "GarfieldMessenger.h"
22#include "G4UIcmdWithADoubleAndUnit.hh"
23#include "G4UIcmdWithAString.hh"
24#include "G4UIcmdWithoutParameter.hh"
25#include "G4UIcommand.hh"
26#include "G4UIdirectory.hh"
27#include "G4UIparameter.hh"
28#include "GarfieldPhysics.h"
29
31 : G4UImessenger(),
32 fExampleDir(0),
33 fGarfieldPhysicsDir(0),
34 fIonizationModelCmd(0),
35 fGarfieldParticleCmd(0)
36{
37 fExampleDir = new G4UIdirectory("/exampleGarfield/");
38 fExampleDir->SetGuidance("Commands specific to this example");
39
40 G4bool broadcast = false;
42 new G4UIdirectory("/exampleGarfield/physics/", broadcast);
43 fGarfieldPhysicsDir->SetGuidance(
44 "Particle and energy ranges for Garfield++ physics model");
45
47 new G4UIcommand("/exampleGarfield/physics/setIonizationModel", this);
48 fIonizationModelCmd->SetGuidance("Select ionization model for Garfield++");
49 fIonizationModelCmd->SetGuidance(
50 " and choose whether to use default particles");
51 fIonizationModelCmd->SetGuidance(" and energy ranges for the chosen model");
52 //
53 G4UIparameter* ionizationModelPrm =
54 new G4UIparameter("ionizationModel", 's', false);
55 ionizationModelPrm->SetGuidance("ionization model (1. PAIPhot, 2. Heed)");
56 ionizationModelPrm->SetGuidance(
57 " 1. Geant4 model, delta electrons transported by Heed");
58 ionizationModelPrm->SetGuidance(" 2. Use directly Heed");
59 fIonizationModelCmd->SetParameter(ionizationModelPrm);
60 //
61 G4UIparameter* useDefaultsPrm = new G4UIparameter("useDefaults", 'b', false);
62 useDefaultsPrm->SetGuidance(
63 "true to use default, false to manually choose particles and energies");
64 fIonizationModelCmd->SetParameter(useDefaultsPrm);
65 //
66 fIonizationModelCmd->AvailableForStates(G4State_PreInit);
67
68 fGarfieldParticleCmd = new G4UIcommand(
69 "/exampleGarfield/physics/setGarfieldParticleTypeAndEnergy", this);
70 fGarfieldParticleCmd->SetGuidance(
71 "Select particle types and energies for Heed model.");
72 fGarfieldParticleCmd->SetGuidance(
73 " For PAI and PAIPhot model choose at which energy electrons are");
74 fGarfieldParticleCmd->SetGuidance(
75 " transported as delta electrons by Heed, and treatment of gammas");
76 //
77 G4UIparameter* particleGarfieldPrm =
78 new G4UIparameter("particleName", 's', false);
79 particleGarfieldPrm->SetGuidance(
80 "Particle name (gamma, e-, e+, mu-, mu+, proton, anti_proton, pi-, pi+, "
81 "kaon, kaon+, alpha, deuteron)");
82 fGarfieldParticleCmd->SetParameter(particleGarfieldPrm);
83 //
84 G4UIparameter* minEnergyGarfieldPrm =
85 new G4UIparameter("minimumEnergyGarfield", 'd', false);
86 minEnergyGarfieldPrm->SetGuidance("minimum energy");
87 minEnergyGarfieldPrm->SetParameterRange("minimumEnergyGarfield>=0");
88 fGarfieldParticleCmd->SetParameter(minEnergyGarfieldPrm);
89 //
90 G4UIparameter* maxEnergyGarfieldPrm =
91 new G4UIparameter("maximumEnergyGarfield", 'd', false);
92 maxEnergyGarfieldPrm->SetGuidance("maximum energy");
93 maxEnergyGarfieldPrm->SetParameterRange("maximumEnergyGarfield>=0");
94 fGarfieldParticleCmd->SetParameter(maxEnergyGarfieldPrm);
95 //
96 G4UIparameter* unitGarfieldPrm = new G4UIparameter("unit", 's', false);
97 unitGarfieldPrm->SetGuidance("unit of energy");
98 G4String unitListGarfield =
99 G4UIcommand::UnitsList(G4UIcommand::CategoryOf("MeV"));
100 unitGarfieldPrm->SetParameterCandidates(unitListGarfield);
101 fGarfieldParticleCmd->SetParameter(unitGarfieldPrm);
102 //
103 fGarfieldParticleCmd->AvailableForStates(G4State_PreInit);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
118void GarfieldMessenger::SetNewValue(G4UIcommand* command, G4String newValue)
119{
120 if (command == fIonizationModelCmd) {
122 G4String modelName;
123 G4bool useDefaults;
124 std::istringstream is(newValue);
125 is >> modelName >> std::boolalpha >> useDefaults;
126 garfieldPhysics->SetIonizationModel(modelName, useDefaults);
127 }
128 else if (command == fGarfieldParticleCmd) {
130 G4String particleName, unit, programName;
131 G4double minEnergy;
132 G4double maxEnergy;
133 std::istringstream is(newValue);
134 is >> particleName >> minEnergy >> maxEnergy >> unit;
135 minEnergy *= G4UIcommand::ValueOf(unit);
136 maxEnergy *= G4UIcommand::ValueOf(unit);
137 garfieldPhysics->AddParticleName(particleName, minEnergy, maxEnergy);
138 }
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the GarfieldMessengerclass.
Definition of the GarfieldPhysics class.
G4UIdirectory * fGarfieldPhysicsDir
G4UIcommand * fGarfieldParticleCmd
G4UIcommand * fIonizationModelCmd
G4UIdirectory * fExampleDir
virtual void SetNewValue(G4UIcommand *, G4String)
static GarfieldPhysics * GetInstance()
void SetIonizationModel(std::string model, bool useDefaults=true)
void AddParticleName(const std::string particleName, double ekin_min_keV, double ekin_max_keV)