#include <BesEmcDetectorMessenger.hh>
Public Member Functions | |
BesEmcDetectorMessenger (BesEmcConstruction *, BesEmcGeometry *) | |
~BesEmcDetectorMessenger () | |
void | SetNewValue (G4UIcommand *, G4String) |
G4String | GetCurrentValue (G4UIcommand *) |
Private Attributes | |
BesEmcConstruction * | BesEmc |
BesEmcGeometry * | fBesEmcGeometry |
G4UIdirectory * | BesdetDir |
G4UIcmdWithAnInteger * | verboseCmd |
G4UIcmdWithAString * | CryMaterCmd |
G4UIcmdWithAString * | CasingMaterCmd |
G4UIcmdWith3VectorAndUnit * | CasingThickCmd |
G4UIcmdWithADoubleAndUnit * | SizeRminCmd |
G4UIcmdWithAnInteger * | NbPhiCmd |
G4UIcmdWithAnInteger * | NbThetaCmd |
G4UIcmdWithAnInteger * | StartIDThetaCmd |
G4UIcmdWithADoubleAndUnit * | LengthCmd |
G4UIcmdWithADoubleAndUnit * | MagFieldCmd |
G4UIcmdWithoutParameter * | UpdateCmd |
Definition at line 27 of file BesEmcDetectorMessenger.hh.
BesEmcDetectorMessenger::BesEmcDetectorMessenger | ( | BesEmcConstruction * | , | |
BesEmcGeometry * | ||||
) |
Definition at line 24 of file BesEmcDetectorMessenger.cc.
References BesdetDir, BesEmc, CasingMaterCmd, CasingThickCmd, CryMaterCmd, fBesEmcGeometry, BesEmcGeometry::GetBSCCryLength(), BesEmcGeometry::GetBSCRmin(), BesEmcGeometry::GetCasingThickness(), BesEmcConstruction::GetMagField(), LengthCmd, MagFieldCmd, NbPhiCmd, NbThetaCmd, SizeRminCmd, StartIDThetaCmd, UpdateCmd, and verboseCmd.
00027 :BesEmc(BesDet) 00028 { 00029 fBesEmcGeometry=besEMCGeometry; 00030 00031 BesdetDir = new G4UIdirectory("/calor/"); 00032 BesdetDir->SetGuidance("Bes calorimeter detector control."); 00033 00034 verboseCmd = new G4UIcmdWithAnInteger("/calor/verbose",this); 00035 verboseCmd->SetGuidance("Set Verbose level of calor management category."); 00036 verboseCmd->SetGuidance(" 0 : Important information"); 00037 verboseCmd->SetGuidance(" 1 : Base information"); 00038 verboseCmd->SetGuidance(" 2 : More..."); 00039 verboseCmd->SetParameterName("level",true); 00040 verboseCmd->SetRange("level>=0"); 00041 verboseCmd->SetDefaultValue(0); 00042 00043 CryMaterCmd = new G4UIcmdWithAString("/calor/setCryMat",this); 00044 CryMaterCmd->SetGuidance("Select Material of the Crystal."); 00045 CryMaterCmd->SetParameterName("choice",true,true); 00046 // CryMaterCmd->AvailableForStates(Idle); 00047 00048 CasingMaterCmd = new G4UIcmdWithAString("/calor/setCasMat",this); 00049 CasingMaterCmd->SetGuidance("Select Material of the Casing."); 00050 CasingMaterCmd->SetParameterName("choice",true,true); 00051 // CasingMaterCmd->AvailableForStates(Idle); 00052 00053 CasingThickCmd = new G4UIcmdWith3VectorAndUnit("/calor/setCasThick",this); 00054 CasingThickCmd->SetGuidance("Set Thickness of the Casing"); 00055 CasingThickCmd->SetParameterName("TyvekThk","AlThk","MylarThk",true); 00056 CasingThickCmd->SetRange("Size>=0."); 00057 CasingThickCmd->SetUnitCategory("Length"); 00058 CasingThickCmd->SetDefaultValue(fBesEmcGeometry->GetCasingThickness()/mm); 00059 CasingThickCmd->SetDefaultUnit("mm"); 00060 // CasingThickCmd->AvailableForStates(Idle); 00061 00062 SizeRminCmd = new G4UIcmdWithADoubleAndUnit("/calor/setSizeRmin",this); 00063 SizeRminCmd->SetGuidance("Set Rmin size of the barrel calorimeter"); 00064 SizeRminCmd->SetParameterName("Size",true); 00065 SizeRminCmd->SetRange("Size>0."); 00066 SizeRminCmd->SetUnitCategory("Length"); 00067 SizeRminCmd->SetDefaultValue(fBesEmcGeometry->GetBSCRmin()/cm); 00068 SizeRminCmd->SetDefaultUnit("cm"); 00069 // SizeRminCmd->AvailableForStates(Idle); 00070 00071 NbPhiCmd = new G4UIcmdWithAnInteger("/calor/setNbPhi",this); 00072 NbPhiCmd->SetGuidance("Set number of crystals at phi direction."); 00073 NbPhiCmd->SetParameterName("NbCrystals",false); 00074 NbPhiCmd->SetRange("NbCrystals>0 && NbCrystals<=150"); 00075 // NbPhiCmd->AvailableForStates(Idle); 00076 00077 NbThetaCmd = new G4UIcmdWithAnInteger("/calor/setNbTheta",this); 00078 NbThetaCmd->SetGuidance("Set number of crystals at theta direction."); 00079 NbThetaCmd->SetParameterName("NbCrystals",false); 00080 NbThetaCmd->SetRange("NbCrystals>0 && NbCrystals<=22"); 00081 // NbThetaCmd->AvailableForStates(Idle); 00082 00083 StartIDThetaCmd = new G4UIcmdWithAnInteger("/calor/setStartTheta",this); 00084 StartIDThetaCmd->SetGuidance("Set ID of starting crystals at theta direction."); 00085 StartIDThetaCmd->SetParameterName("IDCrystal",false); 00086 StartIDThetaCmd->SetRange("IDCrystal>=0 && IDCrystal<22"); 00087 // StartIDThetaCmd->AvailableForStates(Idle); 00088 00089 LengthCmd = new G4UIcmdWithADoubleAndUnit("/calor/setCryLength",this); 00090 LengthCmd->SetGuidance("Set Length of crystals of barrel calorimeter"); 00091 LengthCmd->SetParameterName("Size",true); 00092 LengthCmd->SetRange("Size>0."); 00093 LengthCmd->SetUnitCategory("Length"); 00094 LengthCmd->SetDefaultValue(fBesEmcGeometry->GetBSCCryLength()/cm); 00095 LengthCmd->SetDefaultUnit("cm"); 00096 // LengthCmd->AvailableForStates(Idle); 00097 00098 UpdateCmd = new G4UIcmdWithoutParameter("/calor/update",this); 00099 UpdateCmd->SetGuidance("Update calorimeter geometry."); 00100 UpdateCmd->SetGuidance("This command MUST be applied before \"beamOn\" "); 00101 UpdateCmd->SetGuidance("if you changed geometrical value(s)."); 00102 // UpdateCmd->AvailableForStates(Idle); 00103 00104 MagFieldCmd = new G4UIcmdWithADoubleAndUnit("/calor/setField",this); 00105 MagFieldCmd->SetGuidance("Define magnetic field."); 00106 MagFieldCmd->SetGuidance("Magnetic field will be in Z direction."); 00107 MagFieldCmd->SetParameterName("Bz",true); 00108 MagFieldCmd->SetRange("Bz>=0."); 00109 MagFieldCmd->SetUnitCategory("Magnetic flux density"); 00110 MagFieldCmd->SetDefaultValue(BesEmc->GetMagField()/tesla); 00111 MagFieldCmd->SetDefaultUnit("tesla"); 00112 // MagFieldCmd->AvailableForStates(Idle); 00113 }
BesEmcDetectorMessenger::~BesEmcDetectorMessenger | ( | ) |
Definition at line 117 of file BesEmcDetectorMessenger.cc.
References BesdetDir, CasingMaterCmd, CasingThickCmd, CryMaterCmd, LengthCmd, MagFieldCmd, NbPhiCmd, NbThetaCmd, SizeRminCmd, StartIDThetaCmd, UpdateCmd, and verboseCmd.
00118 { 00119 delete verboseCmd; 00120 delete LengthCmd; 00121 delete NbPhiCmd; delete NbThetaCmd; 00122 delete StartIDThetaCmd; 00123 delete CryMaterCmd; delete CasingMaterCmd; 00124 delete CasingThickCmd; 00125 delete SizeRminCmd; 00126 delete UpdateCmd; 00127 delete MagFieldCmd; 00128 delete BesdetDir; 00129 }
G4String BesEmcDetectorMessenger::GetCurrentValue | ( | G4UIcommand * | ) |
Definition at line 173 of file BesEmcDetectorMessenger.cc.
References BesEmc, CasingMaterCmd, CryMaterCmd, BesEmcConstruction::GetCasingMaterial(), BesEmcConstruction::GetCrystalMaterial(), BesEmcConstruction::GetVerboseLevel(), and verboseCmd.
00174 { 00175 G4String cv; 00176 00177 if( command == verboseCmd ) 00178 { cv = verboseCmd->ConvertToString(BesEmc->GetVerboseLevel()); } 00179 00180 if( command == CryMaterCmd ) 00181 { cv = BesEmc->GetCrystalMaterial()->GetName(); } 00182 00183 if( command == CasingMaterCmd ) 00184 { cv = BesEmc->GetCasingMaterial()->GetName(); } 00185 00186 return cv; 00187 }
void BesEmcDetectorMessenger::SetNewValue | ( | G4UIcommand * | , | |
G4String | ||||
) |
Definition at line 133 of file BesEmcDetectorMessenger.cc.
References BesEmc, CasingMaterCmd, CasingThickCmd, CryMaterCmd, LengthCmd, MagFieldCmd, NbPhiCmd, NbThetaCmd, BesEmcConstruction::SetBSCCrystalLength(), BesEmcConstruction::SetBSCNbPhi(), BesEmcConstruction::SetBSCNbTheta(), BesEmcConstruction::SetBSCRmin(), BesEmcConstruction::SetCasingMaterial(), BesEmcConstruction::SetCasingThickness(), BesEmcConstruction::SetCrystalMaterial(), BesEmcConstruction::SetMagField(), BesEmcConstruction::SetStartIDTheta(), BesEmcConstruction::SetVerboseLevel(), SizeRminCmd, StartIDThetaCmd, UpdateCmd, BesEmcConstruction::UpdateGeometry(), and verboseCmd.
00134 { 00135 if( command == verboseCmd ) 00136 { 00137 BesEmc->SetVerboseLevel(verboseCmd->GetNewIntValue(newValue)); 00138 } 00139 00140 if( command == CryMaterCmd ) 00141 { BesEmc->SetCrystalMaterial(newValue);} 00142 00143 if( command == CasingMaterCmd ) 00144 { BesEmc->SetCasingMaterial(newValue);} 00145 00146 if( command == CasingThickCmd ) 00147 { BesEmc->SetCasingThickness(CasingThickCmd->GetNew3VectorValue(newValue));} 00148 00149 if( command == SizeRminCmd ) 00150 { BesEmc->SetBSCRmin(SizeRminCmd->GetNewDoubleValue(newValue));} 00151 00152 if( command == NbPhiCmd ) 00153 { BesEmc->SetBSCNbPhi(NbPhiCmd->GetNewIntValue(newValue));} 00154 00155 if( command == NbThetaCmd ) 00156 { BesEmc->SetBSCNbTheta(NbThetaCmd->GetNewIntValue(newValue));} 00157 00158 if( command == StartIDThetaCmd ) 00159 { BesEmc->SetStartIDTheta(StartIDThetaCmd->GetNewIntValue(newValue));} 00160 00161 if( command == LengthCmd ) 00162 { BesEmc->SetBSCCrystalLength(LengthCmd->GetNewDoubleValue(newValue));} 00163 00164 if( command == UpdateCmd ) 00165 { BesEmc->UpdateGeometry(); } 00166 00167 if( command == MagFieldCmd ) 00168 { BesEmc->SetMagField(MagFieldCmd->GetNewDoubleValue(newValue));} 00169 }
G4UIdirectory* BesEmcDetectorMessenger::BesdetDir [private] |
Definition at line 39 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), and ~BesEmcDetectorMessenger().
Definition at line 36 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), GetCurrentValue(), and SetNewValue().
G4UIcmdWithAString* BesEmcDetectorMessenger::CasingMaterCmd [private] |
Definition at line 42 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), GetCurrentValue(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWith3VectorAndUnit* BesEmcDetectorMessenger::CasingThickCmd [private] |
Definition at line 43 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithAString* BesEmcDetectorMessenger::CryMaterCmd [private] |
Definition at line 41 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), GetCurrentValue(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithADoubleAndUnit* BesEmcDetectorMessenger::LengthCmd [private] |
Definition at line 48 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithADoubleAndUnit* BesEmcDetectorMessenger::MagFieldCmd [private] |
Definition at line 49 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithAnInteger* BesEmcDetectorMessenger::NbPhiCmd [private] |
Definition at line 45 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithAnInteger* BesEmcDetectorMessenger::NbThetaCmd [private] |
Definition at line 46 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithADoubleAndUnit* BesEmcDetectorMessenger::SizeRminCmd [private] |
Definition at line 44 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithAnInteger* BesEmcDetectorMessenger::StartIDThetaCmd [private] |
Definition at line 47 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithoutParameter* BesEmcDetectorMessenger::UpdateCmd [private] |
Definition at line 50 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), SetNewValue(), and ~BesEmcDetectorMessenger().
G4UIcmdWithAnInteger* BesEmcDetectorMessenger::verboseCmd [private] |
Definition at line 40 of file BesEmcDetectorMessenger.hh.
Referenced by BesEmcDetectorMessenger(), GetCurrentValue(), SetNewValue(), and ~BesEmcDetectorMessenger().