Geant4  10.02.p03
B5DetectorConstruction Class Reference

Detector construction. More...

#include <B5DetectorConstruction.hh>

Inheritance diagram for B5DetectorConstruction:
Collaboration diagram for B5DetectorConstruction:

Public Member Functions

 B5DetectorConstruction ()
 
virtual ~B5DetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
void SetArmAngle (G4double val)
 
G4double GetArmAngle ()
 
void ConstructMaterials ()
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Private Member Functions

void DefineCommands ()
 

Private Attributes

G4GenericMessengerfMessenger
 
G4LogicalVolumefHodoscope1Logical
 
G4LogicalVolumefHodoscope2Logical
 
G4LogicalVolumefWirePlane1Logical
 
G4LogicalVolumefWirePlane2Logical
 
G4LogicalVolumefCellLogical
 
G4LogicalVolumefHadCalScintiLogical
 
G4LogicalVolumefMagneticLogical
 
std::vector< G4VisAttributes * > fVisAttributes
 
G4double fArmAngle
 
G4RotationMatrixfArmRotation
 
G4VPhysicalVolumefSecondArmPhys
 

Static Private Attributes

static G4ThreadLocal B5MagneticFieldfMagneticField = 0
 
static G4ThreadLocal G4FieldManagerfFieldMgr = 0
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 

Detailed Description

Detector construction.

Definition at line 51 of file B5DetectorConstruction.hh.

Constructor & Destructor Documentation

◆ B5DetectorConstruction()

B5DetectorConstruction::B5DetectorConstruction ( )

Definition at line 77 of file B5DetectorConstruction.cc.

79  fMessenger(0),
86 
87 {
90 
91  // define commands for this class
93 }
G4GenericMessenger * fMessenger
G4LogicalVolume * fWirePlane2Logical
CLHEP::HepRotation G4RotationMatrix
G4LogicalVolume * fHodoscope1Logical
G4LogicalVolume * fMagneticLogical
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
G4VPhysicalVolume * fSecondArmPhys
static const double deg
Definition: G4SIunits.hh:151
G4LogicalVolume * fHadCalScintiLogical
G4RotationMatrix * fArmRotation
G4LogicalVolume * fWirePlane1Logical
G4LogicalVolume * fHodoscope2Logical
std::vector< G4VisAttributes * > fVisAttributes
Here is the call graph for this function:

◆ ~B5DetectorConstruction()

B5DetectorConstruction::~B5DetectorConstruction ( )
virtual

Definition at line 97 of file B5DetectorConstruction.cc.

98 {
99  delete fArmRotation;
100  delete fMessenger;
101 
102  for (G4int i=0; i<G4int(fVisAttributes.size()); ++i)
103  {
104  delete fVisAttributes[i];
105  }
106 }
G4GenericMessenger * fMessenger
int G4int
Definition: G4Types.hh:78
G4RotationMatrix * fArmRotation
std::vector< G4VisAttributes * > fVisAttributes

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * B5DetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 110 of file B5DetectorConstruction.cc.

111 {
112  // Construct materials
114  G4Material* air = G4Material::GetMaterial("G4_AIR");
115  //G4Material* argonGas = G4Material::GetMaterial("B5_Ar");
116  G4Material* argonGas = G4Material::GetMaterial("G4_Ar");
117  G4Material* scintillator
118  = G4Material::GetMaterial("G4_PLASTIC_SC_VINYLTOLUENE");
119  G4Material* csI = G4Material::GetMaterial("G4_CESIUM_IODIDE");
120  G4Material* lead = G4Material::GetMaterial("G4_Pb");
121 
122  // Option to switch on/off checking of volumes overlaps
123  //
124  G4bool checkOverlaps = true;
125 
126  // geometries --------------------------------------------------------------
127  // experimental hall (world volume)
128  G4VSolid* worldSolid
129  = new G4Box("worldBox",10.*m,3.*m,10.*m);
130  G4LogicalVolume* worldLogical
131  = new G4LogicalVolume(worldSolid,air,"worldLogical");
132  G4VPhysicalVolume* worldPhysical
133  = new G4PVPlacement(0,G4ThreeVector(),worldLogical,"worldPhysical",0,
134  false,0,checkOverlaps);
135 
136  // Tube with Local Magnetic field
137 
138  G4VSolid* magneticSolid
139  = new G4Tubs("magneticTubs",0.,1.*m,1.*m,0.,360.*deg);
140 
142  = new G4LogicalVolume(magneticSolid, air, "magneticLogical");
143 
144  // placement of Tube
145 
146  G4RotationMatrix* fieldRot = new G4RotationMatrix();
147  fieldRot->rotateX(90.*deg);
149  "magneticPhysical",worldLogical,
150  false,0,checkOverlaps);
151 
152  // set step limit in tube with magnetic field
153  G4UserLimits* userLimits = new G4UserLimits(1*m);
154  fMagneticLogical->SetUserLimits(userLimits);
155 
156  // first arm
157  G4VSolid* firstArmSolid
158  = new G4Box("firstArmBox",1.5*m,1.*m,3.*m);
159  G4LogicalVolume* firstArmLogical
160  = new G4LogicalVolume(firstArmSolid,air,"firstArmLogical");
161  new G4PVPlacement(0,G4ThreeVector(0.,0.,-5.*m),firstArmLogical,
162  "firstArmPhysical",worldLogical,
163  false,0,checkOverlaps);
164 
165  // second arm
166  G4VSolid* secondArmSolid
167  = new G4Box("secondArmBox",2.*m,2.*m,3.5*m);
168  G4LogicalVolume* secondArmLogical
169  = new G4LogicalVolume(secondArmSolid,air,"secondArmLogical");
170  G4double x = -5.*m * std::sin(fArmAngle);
171  G4double z = 5.*m * std::cos(fArmAngle);
173  = new G4PVPlacement(fArmRotation,G4ThreeVector(x,0.,z),secondArmLogical,
174  "fSecondArmPhys",worldLogical,
175  false,0,checkOverlaps);
176 
177  // hodoscopes in first arm
178  G4VSolid* hodoscope1Solid
179  = new G4Box("hodoscope1Box",5.*cm,20.*cm,0.5*cm);
181  = new G4LogicalVolume(hodoscope1Solid,scintillator,"hodoscope1Logical");
182  for (G4int i=0;i<15;i++)
183  {
184  G4double x1 = (i-7)*10.*cm;
186  "hodoscope1Physical",firstArmLogical,
187  false,i,checkOverlaps);
188  }
189 
190  // drift chambers in first arm
191  G4VSolid* chamber1Solid
192  = new G4Box("chamber1Box",1.*m,30.*cm,1.*cm);
193  G4LogicalVolume* chamber1Logical
194  = new G4LogicalVolume(chamber1Solid,argonGas,"chamber1Logical");
195  for (G4int i=0;i<5;i++)
196  {
197  G4double z1 = (i-2)*0.5*m;
198  new G4PVPlacement(0,G4ThreeVector(0.,0.,z1),chamber1Logical,
199  "chamber1Physical",firstArmLogical,
200  false,i,checkOverlaps);
201  }
202 
203  // "virtual" wire plane
204  G4VSolid* wirePlane1Solid
205  = new G4Box("wirePlane1Box",1.*m,30.*cm,0.1*mm);
207  = new G4LogicalVolume(wirePlane1Solid,argonGas,"wirePlane1Logical");
209  "wirePlane1Physical",chamber1Logical,
210  false,0,checkOverlaps);
211 
212  // hodoscopes in second arm
213  G4VSolid* hodoscope2Solid
214  = new G4Box("hodoscope2Box",5.*cm,20.*cm,0.5*cm);
216  = new G4LogicalVolume(hodoscope2Solid,scintillator,"hodoscope2Logical");
217  for (G4int i=0;i<25;i++)
218  {
219  G4double x2 = (i-12)*10.*cm;
221  "hodoscope2Physical",secondArmLogical,
222  false,i,checkOverlaps);
223  }
224 
225  // drift chambers in second arm
226  G4VSolid* chamber2Solid
227  = new G4Box("chamber2Box",1.5*m,30.*cm,1.*cm);
228  G4LogicalVolume* chamber2Logical
229  = new G4LogicalVolume(chamber2Solid,argonGas,"chamber2Logical");
230  for (G4int i=0;i<5;i++)
231  {
232  G4double z2 = (i-2)*0.5*m - 1.5*m;
233  new G4PVPlacement(0,G4ThreeVector(0.,0.,z2),chamber2Logical,
234  "chamber2Physical",secondArmLogical,
235  false,i,checkOverlaps);
236  }
237 
238  // "virtual" wire plane
239  G4VSolid* wirePlane2Solid
240  = new G4Box("wirePlane2Box",1.5*m,30.*cm,0.1*mm);
242  = new G4LogicalVolume(wirePlane2Solid,argonGas,"wirePlane2Logical");
244  "wirePlane2Physical",chamber2Logical,
245  false,0,checkOverlaps);
246 
247  // CsI calorimeter
248  G4VSolid* emCalorimeterSolid
249  = new G4Box("EMcalorimeterBox",1.5*m,30.*cm,15.*cm);
250  G4LogicalVolume* emCalorimeterLogical
251  = new G4LogicalVolume(emCalorimeterSolid,csI,"EMcalorimeterLogical");
252  new G4PVPlacement(0,G4ThreeVector(0.,0.,2.*m),emCalorimeterLogical,
253  "EMcalorimeterPhysical",secondArmLogical,
254  false,0,checkOverlaps);
255 
256  // EMcalorimeter cells
257  G4VSolid* cellSolid
258  = new G4Box("cellBox",7.5*cm,7.5*cm,15.*cm);
260  = new G4LogicalVolume(cellSolid,csI,"cellLogical");
262  new G4PVParameterised("cellPhysical",fCellLogical,emCalorimeterLogical,
263  kXAxis,80,cellParam);
264 
265  // hadron calorimeter
266  G4VSolid* hadCalorimeterSolid
267  = new G4Box("HadCalorimeterBox",1.5*m,30.*cm,50.*cm);
268  G4LogicalVolume* hadCalorimeterLogical
269  = new G4LogicalVolume(hadCalorimeterSolid,lead,"HadCalorimeterLogical");
270  new G4PVPlacement(0,G4ThreeVector(0.,0.,3.*m),hadCalorimeterLogical,
271  "HadCalorimeterPhysical",secondArmLogical,
272  false,0,checkOverlaps);
273 
274  // hadron calorimeter column
275  G4VSolid* HadCalColumnSolid
276  = new G4Box("HadCalColumnBox",15.*cm,30.*cm,50.*cm);
277  G4LogicalVolume* HadCalColumnLogical
278  = new G4LogicalVolume(HadCalColumnSolid,lead,"HadCalColumnLogical");
279  new G4PVReplica("HadCalColumnPhysical",HadCalColumnLogical,
280  hadCalorimeterLogical,kXAxis,10,30.*cm);
281 
282  // hadron calorimeter cell
283  G4VSolid* HadCalCellSolid
284  = new G4Box("HadCalCellBox",15.*cm,15.*cm,50.*cm);
285  G4LogicalVolume* HadCalCellLogical
286  = new G4LogicalVolume(HadCalCellSolid,lead,"HadCalCellLogical");
287  new G4PVReplica("HadCalCellPhysical",HadCalCellLogical,
288  HadCalColumnLogical,kYAxis,2,30.*cm);
289 
290  // hadron calorimeter layers
291  G4VSolid* HadCalLayerSolid
292  = new G4Box("HadCalLayerBox",15.*cm,15.*cm,2.5*cm);
293  G4LogicalVolume* HadCalLayerLogical
294  = new G4LogicalVolume(HadCalLayerSolid,lead,"HadCalLayerLogical");
295  new G4PVReplica("HadCalLayerPhysical",HadCalLayerLogical,
296  HadCalCellLogical,kZAxis,20,5.*cm);
297 
298  // scintillator plates
299  G4VSolid* HadCalScintiSolid
300  = new G4Box("HadCalScintiBox",15.*cm,15.*cm,0.5*cm);
302  = new G4LogicalVolume(HadCalScintiSolid,scintillator,
303  "HadCalScintiLogical");
305  "HadCalScintiPhysical",HadCalLayerLogical,
306  false,0,checkOverlaps);
307 
308  // visualization attributes ------------------------------------------------
309 
310  G4VisAttributes* visAttributes = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
311  visAttributes->SetVisibility(false);
312  worldLogical->SetVisAttributes(visAttributes);
313  fVisAttributes.push_back(visAttributes);
314 
315  visAttributes = new G4VisAttributes(G4Colour(0.9,0.9,0.9)); // LightGray
316  fMagneticLogical->SetVisAttributes(visAttributes);
317  fVisAttributes.push_back(visAttributes);
318 
319  visAttributes = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
320  visAttributes->SetVisibility(false);
321  firstArmLogical->SetVisAttributes(visAttributes);
322  secondArmLogical->SetVisAttributes(visAttributes);
323  fVisAttributes.push_back(visAttributes);
324 
325  visAttributes = new G4VisAttributes(G4Colour(0.8888,0.0,0.0));
326  fHodoscope1Logical->SetVisAttributes(visAttributes);
327  fHodoscope2Logical->SetVisAttributes(visAttributes);
328  fVisAttributes.push_back(visAttributes);
329 
330  visAttributes = new G4VisAttributes(G4Colour(0.0,1.0,0.0));
331  chamber1Logical->SetVisAttributes(visAttributes);
332  chamber2Logical->SetVisAttributes(visAttributes);
333  fVisAttributes.push_back(visAttributes);
334 
335  visAttributes = new G4VisAttributes(G4Colour(0.0,0.8888,0.0));
336  visAttributes->SetVisibility(false);
337  fWirePlane1Logical->SetVisAttributes(visAttributes);
338  fWirePlane2Logical->SetVisAttributes(visAttributes);
339  fVisAttributes.push_back(visAttributes);
340 
341  visAttributes = new G4VisAttributes(G4Colour(0.8888,0.8888,0.0));
342  visAttributes->SetVisibility(false);
343  emCalorimeterLogical->SetVisAttributes(visAttributes);
344  fVisAttributes.push_back(visAttributes);
345 
346  visAttributes = new G4VisAttributes(G4Colour(0.9,0.9,0.0));
347  fCellLogical->SetVisAttributes(visAttributes);
348  fVisAttributes.push_back(visAttributes);
349 
350  visAttributes = new G4VisAttributes(G4Colour(0.0, 0.0, 0.9));
351  hadCalorimeterLogical->SetVisAttributes(visAttributes);
352  fVisAttributes.push_back(visAttributes);
353 
354  visAttributes = new G4VisAttributes(G4Colour(0.0, 0.0, 0.9));
355  visAttributes->SetVisibility(false);
356  HadCalColumnLogical->SetVisAttributes(visAttributes);
357  HadCalCellLogical->SetVisAttributes(visAttributes);
358  HadCalLayerLogical->SetVisAttributes(visAttributes);
359  fHadCalScintiLogical->SetVisAttributes(visAttributes);
360  fVisAttributes.push_back(visAttributes);
361 
362  // return the world physical volume ----------------------------------------
363 
364  return worldPhysical;
365 }
static const double cm
Definition: G4SIunits.hh:118
Double_t x2[nxs]
CLHEP::Hep3Vector G4ThreeVector
G4LogicalVolume * fWirePlane2Logical
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
CLHEP::HepRotation G4RotationMatrix
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Box.hh:64
G4LogicalVolume * fHodoscope1Logical
void SetVisibility(G4bool)
void SetUserLimits(G4UserLimits *pULimits)
Definition: G4Tubs.hh:85
G4LogicalVolume * fMagneticLogical
G4VPhysicalVolume * fSecondArmPhys
int G4int
Definition: G4Types.hh:78
static const double deg
Definition: G4SIunits.hh:151
bool G4bool
Definition: G4Types.hh:79
G4LogicalVolume * fHadCalScintiLogical
EM Calorimeter cell parameterisation.
Double_t x1[nxs]
G4RotationMatrix * fArmRotation
G4LogicalVolume * fWirePlane1Logical
G4LogicalVolume * fHodoscope2Logical
static const double m
Definition: G4SIunits.hh:128
double G4double
Definition: G4Types.hh:76
static const double mm
Definition: G4SIunits.hh:114
void SetVisAttributes(const G4VisAttributes *pVA)
std::vector< G4VisAttributes * > fVisAttributes
Here is the call graph for this function:

◆ ConstructMaterials()

void B5DetectorConstruction::ConstructMaterials ( )

Definition at line 420 of file B5DetectorConstruction.cc.

421 {
422  G4NistManager* nistManager = G4NistManager::Instance();
423 
424  // Air
425  nistManager->FindOrBuildMaterial("G4_AIR");
426 
427  // Argon gas
428  nistManager->FindOrBuildMaterial("G4_Ar");
429  // With a density different from the one defined in NIST
430  // G4double density = 1.782e-03*g/cm3;
431  // nistManager->BuildMaterialWithNewDensity("B5_Ar","G4_Ar",density);
432  // !! cases segmentation fault
433 
434  // Scintillator
435  // (PolyVinylToluene, C_9H_10)
436  nistManager->FindOrBuildMaterial("G4_PLASTIC_SC_VINYLTOLUENE");
437 
438  // CsI
439  nistManager->FindOrBuildMaterial("G4_CESIUM_IODIDE");
440 
441  // Lead
442  nistManager->FindOrBuildMaterial("G4_Pb");
443 
444  // Vacuum "Galactic"
445  // nistManager->FindOrBuildMaterial("G4_Galactic");
446 
447  // Vacuum "Air with low density"
448  // G4Material* air = G4Material::GetMaterial("G4_AIR");
449  // G4double density = 1.0e-5*air->GetDensity();
450  // nistManager
451  // ->BuildMaterialWithNewDensity("Air_lowDensity", "G4_AIR", density);
452 
453  G4cout << G4endl << "The materials defined are : " << G4endl << G4endl;
454  G4cout << *(G4Material::GetMaterialTable()) << G4endl;
455 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConstructSDandField()

void B5DetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 369 of file B5DetectorConstruction.cc.

370 {
371  // sensitive detectors -----------------------------------------------------
373  G4String SDname;
374 
375  G4VSensitiveDetector* hodoscope1
376  = new B5HodoscopeSD(SDname="/hodoscope1");
377  SDman->AddNewDetector(hodoscope1);
379 
380  G4VSensitiveDetector* hodoscope2
381  = new B5HodoscopeSD(SDname="/hodoscope2");
382  SDman->AddNewDetector(hodoscope2);
384 
385  G4VSensitiveDetector* chamber1
386  = new B5DriftChamberSD(SDname="/chamber1");
387  SDman->AddNewDetector(chamber1);
389 
390  G4VSensitiveDetector* chamber2
391  = new B5DriftChamberSD(SDname="/chamber2");
392  SDman->AddNewDetector(chamber2);
394 
395  G4VSensitiveDetector* emCalorimeter
396  = new B5EmCalorimeterSD(SDname="/EMcalorimeter");
397  SDman->AddNewDetector(emCalorimeter);
398  fCellLogical->SetSensitiveDetector(emCalorimeter);
399 
400  G4VSensitiveDetector* hadCalorimeter
401  = new B5HadCalorimeterSD(SDname="/HadCalorimeter");
402  SDman->AddNewDetector(hadCalorimeter);
404 
405  // magnetic field ----------------------------------------------------------
407  fFieldMgr = new G4FieldManager();
410  G4bool forceToAllDaughters = true;
411  fMagneticLogical->SetFieldManager(fFieldMgr, forceToAllDaughters);
412 
413  // Register the field and its manager for deleting
416 }
G4LogicalVolume * fWirePlane2Logical
G4bool SetDetectorField(G4Field *detectorField)
G4LogicalVolume * fHodoscope1Logical
Drift chamber sensitive detector.
G4LogicalVolume * fMagneticLogical
void SetFieldManager(G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
bool G4bool
Definition: G4Types.hh:79
G4LogicalVolume * fHadCalScintiLogical
static G4ThreadLocal B5MagneticField * fMagneticField
Magnetic field.
static G4ThreadLocal G4FieldManager * fFieldMgr
G4LogicalVolume * fWirePlane1Logical
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
G4LogicalVolume * fHodoscope2Logical
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
Hodoscope sensitive detector.
void CreateChordFinder(G4MagneticField *detectorMagField)
EM calorimeter sensitive detector.
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
Hadron calorimeter sensitive detector.
Here is the call graph for this function:

◆ DefineCommands()

void B5DetectorConstruction::DefineCommands ( )
private

Definition at line 480 of file B5DetectorConstruction.cc.

481 {
482  // Define /B5/detector command directory using generic messenger class
483  fMessenger = new G4GenericMessenger(this,
484  "/B5/detector/",
485  "Detector control");
486 
487  // armAngle command
488  G4GenericMessenger::Command& armAngleCmd
489  = fMessenger->DeclareMethodWithUnit("armAngle","deg",
491  "Set rotation angle of the second arm.");
492  armAngleCmd.SetParameterName("angle", true);
493  armAngleCmd.SetRange("angle>=0. && angle<180.");
494  armAngleCmd.SetDefaultValue("30.");
495 }
G4GenericMessenger * fMessenger
This class is generic messenger.
Command & DeclareMethodWithUnit(const G4String &name, const G4String &defaultUnit, const G4AnyMethod &fun, const G4String &doc="")
Command & SetDefaultValue(const G4String &)
Command & SetRange(const G4String &range)
Command & SetParameterName(const G4String &, G4bool, G4bool=false)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetArmAngle()

G4double B5DetectorConstruction::GetArmAngle ( )
inline

Definition at line 61 of file B5DetectorConstruction.hh.

61 { return fArmAngle; }
Here is the call graph for this function:

◆ SetArmAngle()

void B5DetectorConstruction::SetArmAngle ( G4double  val)

Definition at line 459 of file B5DetectorConstruction.cc.

460 {
461  if (!fSecondArmPhys)
462  {
463  G4cerr << "Detector has not yet been constructed." << G4endl;
464  return;
465  }
466 
467  fArmAngle = val;
468  *fArmRotation = G4RotationMatrix(); // make it unit vector
470  G4double x = -5.*m * std::sin(fArmAngle);
471  G4double z = 5.*m * std::cos(fArmAngle);
473 
474  // tell G4RunManager that we change the geometry
476 }
void GeometryHasBeenModified(G4bool prop=true)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
G4VPhysicalVolume * fSecondArmPhys
void SetTranslation(const G4ThreeVector &v)
G4RotationMatrix * fArmRotation
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fArmAngle

G4double B5DetectorConstruction::fArmAngle
private

Definition at line 83 of file B5DetectorConstruction.hh.

◆ fArmRotation

G4RotationMatrix* B5DetectorConstruction::fArmRotation
private

Definition at line 84 of file B5DetectorConstruction.hh.

◆ fCellLogical

G4LogicalVolume* B5DetectorConstruction::fCellLogical
private

Definition at line 77 of file B5DetectorConstruction.hh.

◆ fFieldMgr

G4ThreadLocal G4FieldManager * B5DetectorConstruction::fFieldMgr = 0
staticprivate

Definition at line 71 of file B5DetectorConstruction.hh.

◆ fHadCalScintiLogical

G4LogicalVolume* B5DetectorConstruction::fHadCalScintiLogical
private

Definition at line 78 of file B5DetectorConstruction.hh.

◆ fHodoscope1Logical

G4LogicalVolume* B5DetectorConstruction::fHodoscope1Logical
private

Definition at line 73 of file B5DetectorConstruction.hh.

◆ fHodoscope2Logical

G4LogicalVolume* B5DetectorConstruction::fHodoscope2Logical
private

Definition at line 74 of file B5DetectorConstruction.hh.

◆ fMagneticField

G4ThreadLocal B5MagneticField * B5DetectorConstruction::fMagneticField = 0
staticprivate

Definition at line 70 of file B5DetectorConstruction.hh.

◆ fMagneticLogical

G4LogicalVolume* B5DetectorConstruction::fMagneticLogical
private

Definition at line 79 of file B5DetectorConstruction.hh.

◆ fMessenger

G4GenericMessenger* B5DetectorConstruction::fMessenger
private

Definition at line 68 of file B5DetectorConstruction.hh.

◆ fSecondArmPhys

G4VPhysicalVolume* B5DetectorConstruction::fSecondArmPhys
private

Definition at line 85 of file B5DetectorConstruction.hh.

◆ fVisAttributes

std::vector<G4VisAttributes*> B5DetectorConstruction::fVisAttributes
private

Definition at line 81 of file B5DetectorConstruction.hh.

◆ fWirePlane1Logical

G4LogicalVolume* B5DetectorConstruction::fWirePlane1Logical
private

Definition at line 75 of file B5DetectorConstruction.hh.

◆ fWirePlane2Logical

G4LogicalVolume* B5DetectorConstruction::fWirePlane2Logical
private

Definition at line 76 of file B5DetectorConstruction.hh.


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