Public Member Functions | Private Attributes | List of all members
LBNEDetectorMessenger Class Reference

#include <LBNEDetectorMessenger.hh>

Inheritance diagram for LBNEDetectorMessenger:

Public Member Functions

 LBNEDetectorMessenger (LBNEDetectorConstruction *)
 
 ~LBNEDetectorMessenger ()
 
void SetNewValue (G4UIcommand *, G4String)
 

Private Attributes

LBNEDetectorConstructionLBNEDetector
 
G4UIdirectory * LBNEDir
 
G4UIdirectory * detDir
 
G4UIcmdWithABool * ConstructTarget
 
G4UIcmdWithoutParameter * UpdateCmd
 
G4UIcmdWithoutParameter * ConstructCmd
 
G4UIcmdWithADoubleAndUnit * WaterLayerThickInHorn
 
G4UIcmdWithADoubleAndUnit * SetBeamlineAngle
 
G4UIcmdWithADoubleAndUnit * SetHornCurrent
 
G4UIcmdWithADoubleAndUnit * SetSkinDepthInnerRad
 
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEccentricityIO
 
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEllipticityI
 
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsCurrentEqualizerLongAbsLength
 
std::vector< G4UIcmdWithAString * > SetFileNameFieldMapForCE
 
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerQuadAmpl
 
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerOctAmpl
 
std::vector< G4UIcmdWithADouble * > SetHornsCurrentMultiplier
 
G4UIcmdWithABool * ConstructSimpAbsorber
 
G4UIcmdWithABool * ConstructSculptedAbsorber
 
G4UIcmdWithABool * DisableSpoiler
 
G4UIcmdWithABool * DisableSculptedLayers
 
G4UIcmdWithABool * DisableMask
 
G4UIcmdWithABool * ExpandAlLayers
 
G4UIcmdWithAnInteger * RemoveLayers
 
G4UIcmdWithADoubleAndUnit * DwStrAbsSteelWidth
 
G4UIcmdWithADoubleAndUnit * DwStrAbsConcreteWidth
 
G4UIcmdWithADoubleAndUnit * ZCoordForPerfectFocusing
 

Detailed Description

Definition at line 16 of file LBNEDetectorMessenger.hh.

Constructor & Destructor Documentation

LBNEDetectorMessenger::LBNEDetectorMessenger ( LBNEDetectorConstruction LBNEDet)

Definition at line 2029 of file LBNEDetectorConstruction.cc.

2029  :LBNEDetector(LBNEDet)
2030 {
2031  LBNEDir = new G4UIdirectory("/LBNE/");
2032  LBNEDir->SetGuidance("UI commands for detector geometry");
2033 
2034  detDir = new G4UIdirectory("/LBNE/det/");
2035  detDir->SetGuidance("detector control");
2036 
2037 
2038  ConstructTarget = new G4UIcmdWithABool("/LBNE/det/constructTarget",this);
2039  ConstructTarget->SetParameterName("constructTarget", false);
2040  ConstructTarget->SetGuidance("Target construction on/off");
2041  ConstructTarget->SetParameterName("constructTarget",true);
2042  ConstructTarget->AvailableForStates(G4State_PreInit,G4State_Idle);
2043 
2044  SetBeamlineAngle = new
2045  G4UIcmdWithADoubleAndUnit("/LBNE/det/setBeamlineAngle",this);
2046  SetBeamlineAngle->SetParameterName("angle", false);
2047  SetBeamlineAngle->SetGuidance("Set the angle of the beamline");
2048 
2049  UpdateCmd = new G4UIcmdWithoutParameter("/LBNE/det/update",this);
2050  UpdateCmd->SetGuidance("Update or Construct LBNE geometry. Same difference ");
2051  UpdateCmd->SetGuidance("This command MUST be applied before \"beamOn\" ");
2052  UpdateCmd->SetGuidance("if you changed geometrical value(s).");
2053  UpdateCmd->AvailableForStates(G4State_PreInit);
2054 
2055  ConstructCmd = new G4UIcmdWithoutParameter("/LBNE/det/construct",this);
2056  ConstructCmd->SetGuidance("Construct LBNE geometry. Should be one and only time ");
2057  ConstructCmd->SetGuidance("This command MUST be applied before \"beamOn\" ");
2058  ConstructCmd->SetGuidance("if you changed geometrical value(s).");
2059  ConstructCmd->AvailableForStates(G4State_PreInit);
2060 
2061  //Messenger for simple geometry of absorber
2062  ConstructSimpAbsorber = new G4UIcmdWithABool("/LBNE/det/ConstructSimpAbsorber",this);
2063  ConstructSimpAbsorber->SetParameterName("Simple Absorber",true);
2064  ConstructSimpAbsorber->SetGuidance("Set geometry of Hadron Absorber. If true simulate simple version of absorber");
2065  ConstructSimpAbsorber->SetDefaultValue(false);
2066  ConstructSimpAbsorber->AvailableForStates(G4State_PreInit,G4State_Idle);
2067 
2068  ConstructSculptedAbsorber = new G4UIcmdWithABool("/LBNE/det/ConstructSculptedAbsorber",this);
2069  ConstructSculptedAbsorber->SetParameterName("Aluminum, Sculpted Absorber",true);
2070  ConstructSculptedAbsorber->SetGuidance("Set geometry of Hadron Absorber. If true simulate the 2015 Sculpted absorber");
2071  ConstructSculptedAbsorber->SetDefaultValue(false);
2072  ConstructSculptedAbsorber->AvailableForStates(G4State_PreInit,G4State_Idle);
2073 
2074  DisableSpoiler = new G4UIcmdWithABool("/LBNE/det/DisableSpoiler",this);
2075  DisableSpoiler->SetParameterName("Remove spoiler",true);
2076  DisableSpoiler->SetGuidance("If true simulate the 2015 Sculpted absorber without the spoiler");
2077  DisableSpoiler->SetDefaultValue(false);
2078  DisableSpoiler->AvailableForStates(G4State_PreInit,G4State_Idle);
2079 
2080  DisableSculptedLayers = new G4UIcmdWithABool("/LBNE/det/DisableSculptedLayers",this);
2081  DisableSculptedLayers->SetParameterName("Remove sculpting",true);
2082  DisableSculptedLayers->SetGuidance("If true simulate the 2015 Sculpted absorber without the sculpting on the aluminum layers");
2083  DisableSculptedLayers->SetDefaultValue(false);
2084  DisableSculptedLayers->AvailableForStates(G4State_PreInit,G4State_Idle);
2085 
2086  DisableMask = new G4UIcmdWithABool("/LBNE/det/DisableMask",this);
2087  DisableMask->SetParameterName("Remove sculpting",true);
2088  DisableMask->SetGuidance("If true simulate the 2015 Sculpted absorber without the mask layers between the spoiler and core part of the absorber");
2089  DisableMask->SetDefaultValue(false);
2090  DisableMask->AvailableForStates(G4State_PreInit,G4State_Idle);
2091 
2092  ExpandAlLayers = new G4UIcmdWithABool("/LBNE/det/ExpandAlLayers",this);
2093  ExpandAlLayers->SetParameterName("Expand Sculpted Absorber Al layers to 2.8x2.8 m",true);
2094  ExpandAlLayers->SetGuidance("If true simulate the 2015 Sculpted absorber with a core region to the size of the decay pipe");
2095  ExpandAlLayers->SetDefaultValue(false);
2096  ExpandAlLayers->AvailableForStates(G4State_PreInit,G4State_Idle);
2097 
2098  RemoveLayers = new G4UIcmdWithAnInteger("/LBNE/det/RemoveLayers",this);
2099  RemoveLayers->SetParameterName("Number of layers to remove from the downstream end of the absorber",true);
2100  RemoveLayers->SetGuidance("Simulate the 2015 Sculpted absorber with fewer layers");
2101  RemoveLayers->SetDefaultValue(0);
2102  RemoveLayers->AvailableForStates(G4State_PreInit,G4State_Idle);
2103 
2104  DwStrAbsSteelWidth = new G4UIcmdWithADoubleAndUnit("/LBNE/det/DwStrAbsSteelWidth",this);
2105  DwStrAbsSteelWidth->SetParameterName("Steel Downstr Width",true);
2106  DwStrAbsSteelWidth->SetGuidance("Define the downstream width (m) of the steel volume in the simple geometry of absorber. Default value is ...");
2108  DwStrAbsSteelWidth->SetDefaultUnit("m");
2109  DwStrAbsSteelWidth->AvailableForStates(G4State_PreInit);
2110 
2111  DwStrAbsConcreteWidth = new G4UIcmdWithADoubleAndUnit("/LBNE/det/DwStrAbsConcreteWidth",this);
2112  DwStrAbsConcreteWidth->SetParameterName("Concrete Downstr Width",true);
2113  DwStrAbsConcreteWidth->SetGuidance("Define the downstream width (m) of the concrete volume in the simple geometry of absorber. Default value is ...");
2115  DwStrAbsConcreteWidth ->SetDefaultUnit("m");
2116  DwStrAbsConcreteWidth->AvailableForStates(G4State_PreInit);
2117 
2118  new G4UnitDefinition("kiloampere" , "kA", "Electric current", 1000.*CLHEP::ampere);
2119  SetHornCurrent = new
2120  G4UIcmdWithADoubleAndUnit("/LBNE/det/seHornCurrent",this);
2121  SetHornCurrent->SetParameterName("Horn Current", false);
2122  SetHornCurrent->SetGuidance(" The current for the horn system. (Horn1 and Horn2 are in series.");
2123  SetHornCurrent ->SetDefaultValue(200.0); // CDR, March 2012
2124  SetHornCurrent->SetDefaultUnit("kA");
2125  SetHornCurrent->SetUnitCandidates("kA");
2126  SetHornCurrent->AvailableForStates(G4State_PreInit);
2127  {
2128  SetSkinDepthInnerRad = new G4UIcmdWithADoubleAndUnit("/LBNE/det/SetSkinDepthInnerRad",this);
2129  std::string aGuidance("The skin depth for the inner conductor. \n ");
2130  aGuidance += std::string(" The default assumes infinite. Using Zarko' thesis model, \n");
2131  aGuidance += std::string(" MINOS Document 5694-v1, section 5.2.4, p150 - 156 \n");
2132  SetSkinDepthInnerRad->SetGuidance(aGuidance.c_str());
2133  SetSkinDepthInnerRad->SetParameterName("SetSkinDepthInnerRad",true);
2134  SetSkinDepthInnerRad->SetDefaultValue (1.0e10*CLHEP::mm);
2135  SetSkinDepthInnerRad->AvailableForStates(G4State_PreInit);
2136  }
2137  {
2138  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2139  std::ostringstream aNStrStr; aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "DeltaEccentricityIO";
2140  std::string aNStr(aNStrStr.str());
2141  G4UIcmdWithADoubleAndUnit* aGUICmd = new G4UIcmdWithADoubleAndUnit(aNStr.c_str(),this);
2142  std::ostringstream aNGStrStr; aNGStrStr
2143  << "The Eccentricity between Inner/Outer conductors for Horn number " << (kH+1) << "\n.";
2144  std::string aNGStr(aNGStrStr.str());
2145  aGUICmd->SetGuidance(aNGStr.c_str());
2146  aGUICmd->SetParameterName("SetHornxDeltaEccentricityIO",true);
2147  aGUICmd->SetDefaultValue (0.);
2148  aGUICmd->AvailableForStates(G4State_PreInit);
2149  SetHornsDeltaEccentricityIO.push_back(aGUICmd);
2150  }
2151  }
2152  {
2153  //
2154 
2155  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2156  std::ostringstream aNStrStr; aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "DeltaEllipticityI";
2157  std::string aNStr(aNStrStr.str());
2158  G4UIcmdWithADoubleAndUnit* aGUICmd = new G4UIcmdWithADoubleAndUnit(aNStr.c_str(),this);
2159  std::ostringstream aNGStrStr; aNGStrStr
2160  << "The Ellipticity of the Inner conductor for Horn number " << (kH+1) << "\n.";
2161  std::string aNGStr(aNGStrStr.str());
2162  aGUICmd->SetGuidance(aNGStr.c_str());
2163  aGUICmd->SetParameterName("SetHornxDeltaEllipticityI",true);
2164  aGUICmd->SetDefaultValue (0.);
2165  aGUICmd->AvailableForStates(G4State_PreInit);
2166  SetHornsDeltaEllipticityI.push_back(aGUICmd);
2167  }
2168  }
2169  {
2170  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2171  std::ostringstream aNStrStr;
2172  aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "CurrentMultipiler";
2173  std::string aNStr(aNStrStr.str());
2174  G4UIcmdWithADouble* aGUICmd = new G4UIcmdWithADouble(aNStr.c_str(),this);
2175  std::ostringstream aNGStrStr; aNGStrStr
2176  << "A magic current amplifier to increase/decrease the current for Horn number " << (kH+1) << "\n.";
2177  std::string aNGStr(aNGStrStr.str());
2178  aGUICmd->SetGuidance(aNGStr.c_str());
2179  aGUICmd->SetParameterName("SetHornxCurrentMultiplier",true);
2180  aGUICmd->SetDefaultValue (1.);
2181  aGUICmd->AvailableForStates(G4State_PreInit);
2182  SetHornsCurrentMultiplier.push_back(aGUICmd);
2183  }
2184  }
2185  {
2186  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2187  std::ostringstream aNStrStr;
2188  aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "CurrentEqualizerLongAbsLength";
2189  std::string aNStr(aNStrStr.str());
2190  G4UIcmdWithADoubleAndUnit* aGUICmd = new G4UIcmdWithADoubleAndUnit(aNStr.c_str(),this);
2191  std::ostringstream aNGStrStr; aNGStrStr
2192  << "The charact. length of the decrease of the azimuthal anisotropy \n "
2193  << " due imperfect current equal. section for Horn number " << (kH+1) << "\n.";
2194  std::string aNGStr(aNGStrStr.str());
2195  aGUICmd->SetGuidance(aNGStr.c_str());
2196  aGUICmd->SetParameterName("SetHornxCurrentEqualizerLongAbsLength",true);
2197  aGUICmd->SetDefaultValue (100.);
2198  aGUICmd->AvailableForStates(G4State_PreInit);
2199  SetHornsCurrentEqualizerLongAbsLength.push_back(aGUICmd);
2200  }
2201  }
2202  {
2203  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2204  std::ostringstream aNStrStr;
2205  aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "CurrentEqualizerQuadAmpl";
2206  std::string aNStr(aNStrStr.str());
2207  G4UIcmdWithADouble* aGUICmd = new G4UIcmdWithADouble(aNStr.c_str(),this);
2208  std::ostringstream aNGStrStr; aNGStrStr
2209  << "The amplitude of the azimuthal, quadrupolar, anisotropy \n "
2210  << " due imperfect current equal. section for Horn number " << (kH+1) << "\n.";
2211  std::string aNGStr(aNGStrStr.str());
2212  aGUICmd->SetGuidance(aNGStr.c_str());
2213  aGUICmd->SetParameterName("SetHornxCurrentEqualizerQuadAmpl",true);
2214  aGUICmd->SetDefaultValue (0.);
2215  aGUICmd->AvailableForStates(G4State_PreInit);
2216  SetHornsCurrentEqualizerQuadAmpl.push_back(aGUICmd);
2217  }
2218  }
2219  {
2220  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2221  std::ostringstream aNStrStr;
2222  aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "CurrentEqualizerOctAmpl";
2223  std::string aNStr(aNStrStr.str());
2224  G4UIcmdWithADouble* aGUICmd = new G4UIcmdWithADouble(aNStr.c_str(),this);
2225  std::ostringstream aNGStrStr; aNGStrStr
2226  << "The amplitude of the azimuthal, quadrupolar, anisotropy \n "
2227  << " due imperfect current equal. section for Horn number " << (kH+1) << "\n.";
2228  std::string aNGStr(aNGStrStr.str());
2229  aGUICmd->SetGuidance(aNGStr.c_str());
2230  aGUICmd->SetParameterName("SetHornxCurrentEqualizerOctAmpl",true);
2231  aGUICmd->SetDefaultValue (0.);
2232  aGUICmd->AvailableForStates(G4State_PreInit);
2233  SetHornsCurrentEqualizerOctAmpl.push_back(aGUICmd);
2234  }
2235  }
2236  {
2237  for (size_t kH=0; kH != 3; kH++) { // We will work with only 3 horns. Per executive decision of Laura F.
2238  std::ostringstream aNStrStr;
2239  aNStrStr << "/LBNE/det/SetHorn" << (kH+1) << "CurrentEqualizerMapFileName";
2240  std::string aNStr(aNStrStr.str());
2241  G4UIcmdWithAString* aGUICmd = new G4UIcmdWithAString(aNStr.c_str(),this);
2242  std::ostringstream aNGStrStr; aNGStrStr
2243  << "The file name for the field map created from a priot run, \n "
2244  << " related to imperfect current equal. section for Horn number " << (kH+1) << "\n.";
2245  std::string aNGStr(aNGStrStr.str());
2246  aGUICmd->SetGuidance(aNGStr.c_str());
2247  aGUICmd->SetParameterName("SetHornxCurrentEqualizerFileName",true);
2248  aGUICmd->SetDefaultValue("?");
2249  aGUICmd->AvailableForStates(G4State_PreInit);
2250  SetFileNameFieldMapForCE.push_back(aGUICmd);
2251  }
2252  }
2253  {
2254  ZCoordForPerfectFocusing = new G4UIcmdWithADoubleAndUnit("/LBNE/det/ZCoordForPerfectFocusing",this);
2255  std::string aGuidance("The Z-location at which the perfect focusing process set the track polar angle to zero \n ");
2256  aGuidance += std::string(" If this Z is > 0. , we trigger the instantaneous foucing at this Z \n");
2257  aGuidance += std::string(" else, Z < 0., we trigger the instantaneous foucing at the radius of Target Canister \n");
2258  aGuidance += std::string(" which, for the LBNE 1.2 MW option is 17.6 mm \n");
2259  ZCoordForPerfectFocusing->SetGuidance(aGuidance.c_str());
2260  ZCoordForPerfectFocusing->SetParameterName("ZCoordForPerfectFocusing", true);
2261  ZCoordForPerfectFocusing->SetDefaultUnit ("m");
2262  ZCoordForPerfectFocusing->SetUnitCandidates ("mm cm m");
2263  ZCoordForPerfectFocusing->SetDefaultValue (1.0e103);
2264  ZCoordForPerfectFocusing->AvailableForStates(G4State_PreInit);
2265  }
2266  /*
2267  {
2268  ConstructLBNFNano = new G4UIcmdWithABool("/LBNF/det/ConstructLBNFNano",this);
2269  ConstructLBNFNano->SetParameterName("ConstructLBNFNano", false);
2270  ConstructLBNFNano->SetGuidance("Neutrino Progenitor Spectrometer construction on/off");
2271  ConstructLBNFNano->SetParameterName("ConstructLBNFNano",true);
2272  ConstructLBNFNano->AvailableForStates(G4State_PreInit,G4State_Idle);
2273  }
2274  */
2275 }
G4UIcmdWithADoubleAndUnit * SetHornCurrent
LBNEDetectorConstruction * LBNEDetector
G4UIcmdWithABool * ConstructTarget
std::string string
Definition: nybbler.cc:12
std::vector< G4UIcmdWithAString * > SetFileNameFieldMapForCE
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEccentricityIO
G4UIcmdWithoutParameter * ConstructCmd
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerOctAmpl
G4UIcmdWithADoubleAndUnit * SetSkinDepthInnerRad
G4UIcmdWithADoubleAndUnit * ZCoordForPerfectFocusing
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEllipticityI
std::vector< G4UIcmdWithADouble * > SetHornsCurrentMultiplier
G4UIcmdWithADoubleAndUnit * SetBeamlineAngle
G4UIcmdWithADoubleAndUnit * DwStrAbsConcreteWidth
G4UIcmdWithABool * ExpandAlLayers
G4UIcmdWithAnInteger * RemoveLayers
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsCurrentEqualizerLongAbsLength
G4UIcmdWithABool * DisableSpoiler
G4UIcmdWithoutParameter * UpdateCmd
G4UIcmdWithABool * ConstructSculptedAbsorber
G4UIcmdWithABool * ConstructSimpAbsorber
G4UIcmdWithABool * DisableSculptedLayers
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerQuadAmpl
G4UIcmdWithABool * DisableMask
G4UIcmdWithADoubleAndUnit * DwStrAbsSteelWidth
LBNEDetectorMessenger::~LBNEDetectorMessenger ( )

Definition at line 2277 of file LBNEDetectorConstruction.cc.

2278 {
2279 
2280  delete detDir;
2281  delete LBNEDir;
2282  delete ConstructSimpAbsorber;
2284  delete DwStrAbsConcreteWidth;
2285  delete DwStrAbsSteelWidth;
2286  delete ConstructTarget;
2287  delete SetBeamlineAngle;
2288  delete UpdateCmd;
2289  delete ConstructCmd;
2290  delete SetHornCurrent;
2291  delete SetSkinDepthInnerRad;
2292  for (size_t kH=0; kH != 3; kH++) {
2293  delete SetHornsDeltaEccentricityIO[kH];
2294  delete SetHornsDeltaEllipticityI[kH];
2297  delete SetFileNameFieldMapForCE[kH];
2299  delete SetHornsCurrentMultiplier[kH];
2300  }
2301  delete ZCoordForPerfectFocusing;
2302  delete DisableSpoiler;
2303  delete DisableSculptedLayers;
2304  delete DisableMask;
2305  delete ExpandAlLayers;
2306  delete RemoveLayers;
2307 }
G4UIcmdWithADoubleAndUnit * SetHornCurrent
G4UIcmdWithABool * ConstructTarget
std::vector< G4UIcmdWithAString * > SetFileNameFieldMapForCE
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEccentricityIO
G4UIcmdWithoutParameter * ConstructCmd
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerOctAmpl
G4UIcmdWithADoubleAndUnit * SetSkinDepthInnerRad
G4UIcmdWithADoubleAndUnit * ZCoordForPerfectFocusing
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEllipticityI
std::vector< G4UIcmdWithADouble * > SetHornsCurrentMultiplier
G4UIcmdWithADoubleAndUnit * SetBeamlineAngle
G4UIcmdWithADoubleAndUnit * DwStrAbsConcreteWidth
G4UIcmdWithABool * ExpandAlLayers
G4UIcmdWithAnInteger * RemoveLayers
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsCurrentEqualizerLongAbsLength
G4UIcmdWithABool * DisableSpoiler
G4UIcmdWithoutParameter * UpdateCmd
G4UIcmdWithABool * ConstructSculptedAbsorber
G4UIcmdWithABool * ConstructSimpAbsorber
G4UIcmdWithABool * DisableSculptedLayers
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerQuadAmpl
G4UIcmdWithABool * DisableMask
G4UIcmdWithADoubleAndUnit * DwStrAbsSteelWidth

Member Function Documentation

void LBNEDetectorMessenger::SetNewValue ( G4UIcommand *  command,
G4String  newValue 
)

Definition at line 2310 of file LBNEDetectorConstruction.cc.

2311 {
2312 
2313  //Turn on and off simple geometry
2314  //Note: ConstructSimpAbsorber defined in LBNEDetectorMessemger.hh
2315  // SetConstructSimpAbsorber and fConstructSimpAbsorber defined in LBNEDetectorConstruction.hh
2316  if ( command == ConstructSimpAbsorber )
2317  {
2318  LBNEDetector->SetConstructSimpAbsorber(ConstructSimpAbsorber->GetNewBoolValue(newValue));
2319  }
2320 
2322  {
2324  }
2325 
2326  if ( command == DisableSpoiler )
2327  {
2328  LBNEDetector->SetDisableSpoiler(DisableSpoiler->GetNewBoolValue(newValue));
2329  }
2330 
2331  if ( command == DisableSculptedLayers )
2332  {
2333  LBNEDetector->SetDisableSculptedLayers(DisableSculptedLayers->GetNewBoolValue(newValue));
2334  }
2335 
2336  if ( command == DisableMask )
2337  {
2338  LBNEDetector->SetDisableMask(DisableMask->GetNewBoolValue(newValue));
2339  }
2340 
2341  if ( command == ExpandAlLayers )
2342  {
2343  LBNEDetector->SetExpandAlLayers(ExpandAlLayers->GetNewBoolValue(newValue));
2344  }
2345 
2346  if ( command == RemoveLayers )
2347  {
2348  LBNEDetector->SetRemoveLayers(RemoveLayers->GetNewIntValue(newValue));
2349  }
2350 
2351  if ( command == DwStrAbsSteelWidth )
2352  {
2353  G4UIcmdWithADoubleAndUnit* cmdSteelDwStr = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2354  LBNEDetector->SetDwStrAbsSteelWidth(cmdSteelDwStr->GetNewDoubleValue(newValue));
2355  //G4double SteelWidth = LBNEDetector->GetDwStrAbsSteelWidth();
2356  }
2357 
2358 
2360  {
2361  G4UIcmdWithADoubleAndUnit* cmdConcreteDwStr = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2362  LBNEDetector->SetDwStrAbsConcreteWidth(cmdConcreteDwStr->GetNewDoubleValue(newValue));
2363  // G4double ConcWidth = LBNEDetector->GetDwStrAbsConcreteWidth();
2364  }
2365 
2366  if ( command == UpdateCmd )
2367  {
2369  return;
2370  }
2371  if ( command == ConstructCmd )
2372  {
2374  return;
2375  }
2376  if (command == SetHornCurrent) {
2377  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2378  LBNEDetector->SetHornCurrent(cmdWD->GetNewDoubleValue(newValue));
2379  }
2380  if (command == SetSkinDepthInnerRad) {
2381  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2382  LBNEDetector->SetSkinDepthInnerRad(cmdWD->GetNewDoubleValue(newValue));
2383  }
2384  for (size_t kH=0; kH != 3; kH++) {
2385  if (command == SetHornsDeltaEccentricityIO[kH]) {
2386  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2387  LBNEDetector->SetDeltaEccentricityIO(kH, cmdWD->GetNewDoubleValue(newValue));
2388  }
2389  if (command == SetHornsDeltaEllipticityI[kH]) {
2390  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2391  LBNEDetector->SetDeltaEllipticityI(kH, cmdWD->GetNewDoubleValue(newValue));
2392  }
2394  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2395  LBNEDetector->SetCurrentEqualizerLongAbsLength(kH, cmdWD->GetNewDoubleValue(newValue));
2396  }
2398  G4UIcmdWithADouble* cmdWD = dynamic_cast<G4UIcmdWithADouble*> (command);
2399  LBNEDetector->SetCurrentEqualizerQuadAmpl(kH, cmdWD->GetNewDoubleValue(newValue));
2400  }
2401  if (command == SetFileNameFieldMapForCE[kH]) {
2402  LBNEDetector->SetFileNameFieldMapForCE(kH, newValue);
2403  }
2405  G4UIcmdWithADouble* cmdWD = dynamic_cast<G4UIcmdWithADouble*> (command);
2406  LBNEDetector->SetCurrentEqualizerOctAmpl(kH, cmdWD->GetNewDoubleValue(newValue));
2407  }
2408  if (command == SetHornsCurrentMultiplier[kH]) {
2409  G4UIcmdWithADouble* cmdWD = dynamic_cast<G4UIcmdWithADouble*> (command);
2410  LBNEDetector->SetCurrentMultipiler(kH, cmdWD->GetNewDoubleValue(newValue));
2411  }
2412  }
2414  G4UIcmdWithADoubleAndUnit* cmdWD = dynamic_cast<G4UIcmdWithADoubleAndUnit*> (command);
2415  double val = cmdWD->GetNewDoubleValue(newValue);
2417  if (val < 0.) {
2418  std::cout << " You have chosen a near perfect focus exiting the target Helium containment vessel " << std::endl;
2419  std::cout << " This only make sens if the real horn current is vanishlingly small.. " << std::endl;
2420  LBNEDetector->SetHornCurrent(2.0e-11);
2421  std::cout << " ...Done.. Also, Horn's conductors are made of air " << std::endl;
2423  G4String matAir("Air");
2424  volP->SetHorn1InnerConductorMaterial(matAir);
2425  volP->SetHorn1AllConductorMaterial(matAir);
2426  volP->SetHorn2InnerConductorMaterial(matAir);
2427  volP->SetHorn2AllConductorMaterial(matAir);
2428  } else if (val > 22499.*CLHEP::mm) {
2429  std::cout << " You have chosen a near perfect focus entering the decay channel at Z = " << val << std::endl;
2430  }
2431  }
2432  /*
2433  ** not yet deployed..
2434  **
2435  if (command == ConstructLBNFNano) {
2436  G4UIcmdWithABool* cmdWD = dynamic_cast<G4UIcmdWithABool*> (command);
2437  LBNEDetector->SetConstructLBNFNano(cmdWD->GetNewBoolValue(newValue));
2438  }
2439  */
2440 }
G4UIcmdWithADoubleAndUnit * SetHornCurrent
void SetConstructSimpAbsorber(G4bool aBool)
LBNEDetectorConstruction * LBNEDetector
std::vector< G4UIcmdWithAString * > SetFileNameFieldMapForCE
static LBNEVolumePlacements * Instance()
int command
void SetCurrentEqualizerOctAmpl(size_t kH, double v)
void SetCurrentEqualizerQuadAmpl(size_t kH, double v)
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEccentricityIO
G4UIcmdWithoutParameter * ConstructCmd
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerOctAmpl
G4UIcmdWithADoubleAndUnit * SetSkinDepthInnerRad
G4UIcmdWithADoubleAndUnit * ZCoordForPerfectFocusing
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsDeltaEllipticityI
std::vector< G4UIcmdWithADouble * > SetHornsCurrentMultiplier
void SetDeltaEccentricityIO(size_t kH, G4double val)
void SetHorn2AllConductorMaterial(G4String &name)
G4UIcmdWithADoubleAndUnit * DwStrAbsConcreteWidth
G4UIcmdWithABool * ExpandAlLayers
void SetHorn2InnerConductorMaterial(G4String &name)
G4UIcmdWithAnInteger * RemoveLayers
void SetDisableSculptedLayers(G4bool aBool)
void SetConstructSculptedAbsorber(G4bool aBool)
std::vector< G4UIcmdWithADoubleAndUnit * > SetHornsCurrentEqualizerLongAbsLength
void SetHorn1InnerConductorMaterial(G4String &name)
G4UIcmdWithABool * DisableSpoiler
G4UIcmdWithoutParameter * UpdateCmd
void SetCurrentEqualizerLongAbsLength(size_t kH, double v)
G4UIcmdWithABool * ConstructSculptedAbsorber
void SetFileNameFieldMapForCE(size_t kH, std::string s)
void SetDeltaEllipticityI(size_t kH, G4double val)
G4UIcmdWithABool * ConstructSimpAbsorber
G4UIcmdWithABool * DisableSculptedLayers
std::vector< G4UIcmdWithADouble * > SetHornsCurrentEqualizerQuadAmpl
G4UIcmdWithABool * DisableMask
G4UIcmdWithADoubleAndUnit * DwStrAbsSteelWidth
void SetCurrentMultipiler(size_t kH, double v)
QTextStream & endl(QTextStream &s)
void SetHorn1AllConductorMaterial(G4String &name)

Member Data Documentation

G4UIcmdWithoutParameter* LBNEDetectorMessenger::ConstructCmd
private

Definition at line 35 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::ConstructSculptedAbsorber
private

Definition at line 52 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::ConstructSimpAbsorber
private

Definition at line 51 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::ConstructTarget
private

Definition at line 31 of file LBNEDetectorMessenger.hh.

G4UIdirectory* LBNEDetectorMessenger::detDir
private

Definition at line 28 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::DisableMask
private

Definition at line 55 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::DisableSculptedLayers
private

Definition at line 54 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::DisableSpoiler
private

Definition at line 53 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::DwStrAbsConcreteWidth
private

Definition at line 60 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::DwStrAbsSteelWidth
private

Definition at line 59 of file LBNEDetectorMessenger.hh.

G4UIcmdWithABool* LBNEDetectorMessenger::ExpandAlLayers
private

Definition at line 56 of file LBNEDetectorMessenger.hh.

LBNEDetectorConstruction* LBNEDetectorMessenger::LBNEDetector
private

Definition at line 25 of file LBNEDetectorMessenger.hh.

G4UIdirectory* LBNEDetectorMessenger::LBNEDir
private

Definition at line 27 of file LBNEDetectorMessenger.hh.

G4UIcmdWithAnInteger* LBNEDetectorMessenger::RemoveLayers
private

Definition at line 57 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::SetBeamlineAngle
private

Definition at line 40 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithAString*> LBNEDetectorMessenger::SetFileNameFieldMapForCE
private

Definition at line 46 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::SetHornCurrent
private

Definition at line 41 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADoubleAndUnit*> LBNEDetectorMessenger::SetHornsCurrentEqualizerLongAbsLength
private

Definition at line 45 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADouble*> LBNEDetectorMessenger::SetHornsCurrentEqualizerOctAmpl
private

Definition at line 48 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADouble*> LBNEDetectorMessenger::SetHornsCurrentEqualizerQuadAmpl
private

Definition at line 47 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADouble*> LBNEDetectorMessenger::SetHornsCurrentMultiplier
private

Definition at line 49 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADoubleAndUnit*> LBNEDetectorMessenger::SetHornsDeltaEccentricityIO
private

Definition at line 43 of file LBNEDetectorMessenger.hh.

std::vector<G4UIcmdWithADoubleAndUnit*> LBNEDetectorMessenger::SetHornsDeltaEllipticityI
private

Definition at line 44 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::SetSkinDepthInnerRad
private

Definition at line 42 of file LBNEDetectorMessenger.hh.

G4UIcmdWithoutParameter* LBNEDetectorMessenger::UpdateCmd
private

Definition at line 34 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::WaterLayerThickInHorn
private

Definition at line 37 of file LBNEDetectorMessenger.hh.

G4UIcmdWithADoubleAndUnit* LBNEDetectorMessenger::ZCoordForPerfectFocusing
private

Definition at line 63 of file LBNEDetectorMessenger.hh.


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