File indexing completed on 2025-02-23 09:22:35
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026 #include "Par04DetectorMessenger.hh"
0027
0028 #include "Par04DetectorConstruction.hh" // for Par04DetectorConstruction
0029
0030 #include "G4UIcmdWithADoubleAndUnit.hh" // for G4UIcmdWithADoubleAndUnit
0031 #include "G4UIcmdWithAnInteger.hh" // for G4UIcmdWithAnInteger
0032 #include "G4UIcmdWithoutParameter.hh" // for G4UIcmdWithoutParameter
0033 #include "G4UIdirectory.hh" // for G4UIdirectory
0034
0035 #include <CLHEP/Units/SystemOfUnits.h> // for pi
0036 #include <G4ApplicationState.hh> // for G4State_PreInit, G4State_Idle
0037 #include <G4ThreeVector.hh> // for G4ThreeVector
0038 #include <G4Types.hh> // for G4bool, G4double, G4int
0039 #include <G4UIcommand.hh> // for G4UIcommand
0040 #include <G4UImessenger.hh> // for G4UImessenger
0041 #include <G4UIparameter.hh> // for G4UIparameter
0042 #include <istream> // for basic_istream, basic_istream...
0043 #include <string> // for operator>>
0044
0045
0046
0047 Par04DetectorMessenger::Par04DetectorMessenger(Par04DetectorConstruction* aDetector)
0048 : G4UImessenger(), fDetector(aDetector)
0049 {
0050 fExampleDir = new G4UIdirectory("/Par04/");
0051 fExampleDir->SetGuidance("UI commands specific to this example");
0052
0053 fDetectorDir = new G4UIdirectory("/Par04/detector/");
0054 fDetectorDir->SetGuidance("Detector construction UI commands");
0055
0056 fPrintCmd = new G4UIcmdWithoutParameter("/Par04/detector/print", this);
0057 fPrintCmd->SetGuidance("Print current settings.");
0058
0059 fDetectorInnerRadiusCmd =
0060 new G4UIcmdWithADoubleAndUnit("/Par04/detector/setDetectorInnerRadius", this);
0061 fDetectorInnerRadiusCmd->SetGuidance("Set cylindrical detector inner radius");
0062 fDetectorInnerRadiusCmd->SetParameterName("Size", false);
0063 fDetectorInnerRadiusCmd->SetRange("Size>0.");
0064 fDetectorInnerRadiusCmd->SetUnitCategory("Length");
0065 fDetectorInnerRadiusCmd->AvailableForStates(G4State_PreInit);
0066 fDetectorInnerRadiusCmd->SetToBeBroadcasted(false);
0067
0068 fDetectorLengthCmd = new G4UIcmdWithADoubleAndUnit("/Par04/detector/setDetectorLength", this);
0069 fDetectorLengthCmd->SetGuidance("Set length of the detector (cylinder length)");
0070 fDetectorLengthCmd->SetParameterName("Size", false);
0071 fDetectorLengthCmd->SetRange("Size>0.");
0072 fDetectorLengthCmd->SetUnitCategory("Length");
0073 fDetectorLengthCmd->AvailableForStates(G4State_PreInit);
0074 fDetectorLengthCmd->SetToBeBroadcasted(false);
0075
0076 fNbLayersCmd = new G4UIcmdWithAnInteger("/Par04/detector/setNbOfLayers", this);
0077 fNbLayersCmd->SetGuidance("Set number of layers.");
0078 fNbLayersCmd->SetParameterName("NbLayers", false);
0079 fNbLayersCmd->SetRange("NbLayers>0");
0080 fNbLayersCmd->AvailableForStates(G4State_PreInit);
0081 fNbLayersCmd->SetToBeBroadcasted(false);
0082
0083 fAbsorCmd = new G4UIcommand("/Par04/detector/setAbsorber", this);
0084 fAbsorCmd->SetGuidance("Set the absorber id, the material, the thickness.");
0085 fAbsorCmd->SetGuidance(" absorber number : from 0 to 1");
0086 fAbsorCmd->SetGuidance(" material name");
0087 fAbsorCmd->SetGuidance(" thickness (with unit) : t>0");
0088 fAbsorCmd->SetGuidance(" if sensitive : true/false.");
0089 auto absNbPrm = new G4UIparameter("AbsorNb", 'i', false);
0090 absNbPrm->SetGuidance("absor number : from 0 to 1");
0091 absNbPrm->SetParameterRange("AbsorNb>-1&AbsoNb<2");
0092 fAbsorCmd->SetParameter(absNbPrm);
0093 auto matPrm = new G4UIparameter("material", 's', false);
0094 matPrm->SetGuidance("material name");
0095 fAbsorCmd->SetParameter(matPrm);
0096 auto thickPrm = new G4UIparameter("thickness", 'd', false);
0097 thickPrm->SetGuidance("thickness of absorber");
0098 thickPrm->SetParameterRange("thickness>0.");
0099 fAbsorCmd->SetParameter(thickPrm);
0100 auto unitPrm = new G4UIparameter("unit", 's', false);
0101 unitPrm->SetGuidance("unit of thickness");
0102 G4String unitList = G4UIcommand::UnitsList(G4UIcommand::CategoryOf("mm"));
0103 unitPrm->SetParameterCandidates(unitList);
0104 fAbsorCmd->SetParameter(unitPrm);
0105 auto sensitivePrm = new G4UIparameter("sensitive", 'b', false);
0106 sensitivePrm->SetGuidance("if absorber is sensitive (registers energy deposits)");
0107 fAbsorCmd->SetParameter(sensitivePrm);
0108
0109 fAbsorCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
0110 fAbsorCmd->SetToBeBroadcasted(false);
0111
0112 fMeshDir = new G4UIdirectory("/Par04/mesh/");
0113 fMeshDir->SetGuidance("Mesh UI commands");
0114
0115 fMeshNbRhoCellsCmd = new G4UIcmdWithAnInteger("/Par04/mesh/setNbOfRhoCells", this);
0116 fMeshNbRhoCellsCmd->SetGuidance("Set number of rho cells in the cylindrical mesh readout.");
0117 fMeshNbRhoCellsCmd->SetParameterName("NbRhoCells", false);
0118 fMeshNbRhoCellsCmd->SetRange("NbRhoCells>0");
0119 fMeshNbRhoCellsCmd->AvailableForStates(G4State_PreInit);
0120 fMeshNbRhoCellsCmd->SetToBeBroadcasted(false);
0121
0122 fMeshNbPhiCellsCmd = new G4UIcmdWithAnInteger("/Par04/mesh/setNbOfPhiCells", this);
0123 fMeshNbPhiCellsCmd->SetGuidance("Set number of phi cells in the cylindrical mesh readout.");
0124 fMeshNbPhiCellsCmd->SetParameterName("NbPhiCells", false);
0125 fMeshNbPhiCellsCmd->SetRange("NbPhiCells>0");
0126 fMeshNbPhiCellsCmd->AvailableForStates(G4State_PreInit);
0127 fMeshNbPhiCellsCmd->SetToBeBroadcasted(false);
0128
0129 fMeshNbZCellsCmd = new G4UIcmdWithAnInteger("/Par04/mesh/setNbOfZCells", this);
0130 fMeshNbZCellsCmd->SetGuidance("Set number of z cells in the cylindrical mesh readout.");
0131 fMeshNbZCellsCmd->SetParameterName("NbZCells", false);
0132 fMeshNbZCellsCmd->SetRange("NbZCells>0");
0133 fMeshNbZCellsCmd->AvailableForStates(G4State_PreInit);
0134 fMeshNbZCellsCmd->SetToBeBroadcasted(false);
0135
0136 fMeshSizeRhoCellsCmd = new G4UIcmdWithADoubleAndUnit("/Par04/mesh/setSizeOfRhoCells", this);
0137 fMeshSizeRhoCellsCmd->SetGuidance("Set size of rho cells in the cylindrical readout mesh");
0138 fMeshSizeRhoCellsCmd->SetParameterName("Size", false);
0139 fMeshSizeRhoCellsCmd->SetRange("Size>0.");
0140 fMeshSizeRhoCellsCmd->SetUnitCategory("Length");
0141 fMeshSizeRhoCellsCmd->AvailableForStates(G4State_PreInit);
0142 fMeshSizeRhoCellsCmd->SetToBeBroadcasted(false);
0143
0144 fMeshSizeZCellsCmd = new G4UIcmdWithADoubleAndUnit("/Par04/mesh/setSizeOfZCells", this);
0145 fMeshSizeZCellsCmd->SetGuidance("Set size of z cells in the cylindrical readout mesh");
0146 fMeshSizeZCellsCmd->SetParameterName("Size", false);
0147 fMeshSizeZCellsCmd->SetRange("Size>0.");
0148 fMeshSizeZCellsCmd->SetUnitCategory("Length");
0149 fMeshSizeZCellsCmd->AvailableForStates(G4State_PreInit);
0150 fMeshSizeZCellsCmd->SetToBeBroadcasted(false);
0151 }
0152
0153
0154
0155 Par04DetectorMessenger::~Par04DetectorMessenger()
0156 {
0157 delete fPrintCmd;
0158 delete fDetectorInnerRadiusCmd;
0159 delete fDetectorLengthCmd;
0160 delete fNbLayersCmd;
0161 delete fAbsorCmd;
0162 delete fDetectorDir;
0163 delete fMeshNbRhoCellsCmd;
0164 delete fMeshNbPhiCellsCmd;
0165 delete fMeshNbZCellsCmd;
0166 delete fMeshSizeRhoCellsCmd;
0167 delete fMeshSizeZCellsCmd;
0168 delete fMeshDir;
0169 delete fExampleDir;
0170 }
0171
0172
0173
0174 void Par04DetectorMessenger::SetNewValue(G4UIcommand* aCommand, G4String aNewValue)
0175 {
0176 if (aCommand == fPrintCmd) {
0177 fDetector->Print();
0178 }
0179 else if (aCommand == fDetectorInnerRadiusCmd) {
0180 fDetector->SetInnerRadius(fDetectorInnerRadiusCmd->GetNewDoubleValue(aNewValue));
0181 }
0182 else if (aCommand == fDetectorLengthCmd) {
0183 fDetector->SetLength(fDetectorInnerRadiusCmd->GetNewDoubleValue(aNewValue));
0184 }
0185 else if (aCommand == fNbLayersCmd) {
0186 fDetector->SetNbOfLayers(fNbLayersCmd->GetNewIntValue(aNewValue));
0187 }
0188 else if (aCommand == fAbsorCmd) {
0189 G4int num;
0190 G4double thick;
0191 G4String unt, mat;
0192 G4bool sensitive;
0193 std::istringstream is(aNewValue);
0194 is >> num >> mat >> thick >> unt >> std::boolalpha >> sensitive;
0195 G4String material = mat;
0196 thick *= G4UIcommand::ValueOf(unt);
0197 fDetector->SetAbsorberMaterial(num, material);
0198 fDetector->SetAbsorberThickness(num, thick);
0199 fDetector->SetAbsorberSensitivity(num, sensitive);
0200 }
0201 else if (aCommand == fMeshNbRhoCellsCmd) {
0202 fDetector->SetMeshNbOfCells(0, fMeshNbRhoCellsCmd->GetNewIntValue(aNewValue));
0203 }
0204 else if (aCommand == fMeshNbPhiCellsCmd) {
0205 fDetector->SetMeshNbOfCells(1, fMeshNbPhiCellsCmd->GetNewIntValue(aNewValue));
0206 fDetector->SetMeshSizeOfCells(1,
0207 2. * CLHEP::pi / fMeshNbPhiCellsCmd->GetNewIntValue(aNewValue));
0208 }
0209 else if (aCommand == fMeshNbZCellsCmd) {
0210 fDetector->SetMeshNbOfCells(2, fMeshNbZCellsCmd->GetNewIntValue(aNewValue));
0211 }
0212 else if (aCommand == fMeshSizeRhoCellsCmd) {
0213 fDetector->SetMeshSizeOfCells(0, fMeshSizeRhoCellsCmd->GetNewDoubleValue(aNewValue));
0214 }
0215 else if (aCommand == fMeshSizeZCellsCmd) {
0216 fDetector->SetMeshSizeOfCells(2, fMeshSizeZCellsCmd->GetNewDoubleValue(aNewValue));
0217 }
0218 }
0219
0220
0221
0222 G4String Par04DetectorMessenger::GetCurrentValue(G4UIcommand* aCommand)
0223 {
0224 G4String cv;
0225
0226 if (aCommand == fDetectorInnerRadiusCmd) {
0227 cv = fDetectorInnerRadiusCmd->ConvertToString(fDetector->GetInnerRadius(), "mm");
0228 }
0229 else if (aCommand == fDetectorLengthCmd) {
0230 cv = fDetectorLengthCmd->ConvertToString(fDetector->GetLength(), "mm");
0231 }
0232 else if (aCommand == fNbLayersCmd) {
0233 cv = fNbLayersCmd->ConvertToString(fDetector->GetNbOfLayers());
0234 }
0235 else if (aCommand == fMeshNbRhoCellsCmd) {
0236 cv = fMeshNbRhoCellsCmd->ConvertToString(fDetector->GetMeshNbOfCells()[0]);
0237 }
0238 else if (aCommand == fMeshNbPhiCellsCmd) {
0239 cv = fMeshNbPhiCellsCmd->ConvertToString(fDetector->GetMeshNbOfCells()[1]);
0240 }
0241 else if (aCommand == fMeshNbZCellsCmd) {
0242 cv = fMeshNbZCellsCmd->ConvertToString(fDetector->GetMeshNbOfCells()[2]);
0243 }
0244 else if (aCommand == fMeshSizeRhoCellsCmd) {
0245 cv = fMeshSizeRhoCellsCmd->ConvertToString(fDetector->GetMeshSizeOfCells()[0]);
0246 }
0247 else if (aCommand == fMeshSizeZCellsCmd) {
0248 cv = fMeshSizeZCellsCmd->ConvertToString(fDetector->GetMeshSizeOfCells()[2]);
0249 }
0250 return cv;
0251 }