Geant4  10.02.p03
B4dDetectorConstruction Class Reference

#include <B4dDetectorConstruction.hh>

Inheritance diagram for B4dDetectorConstruction:
Collaboration diagram for B4dDetectorConstruction:

Public Member Functions

 B4dDetectorConstruction ()
 
virtual ~B4dDetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
- 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 DefineMaterials ()
 
G4VPhysicalVolumeDefineVolumes ()
 

Private Attributes

G4bool fCheckOverlaps
 

Static Private Attributes

static G4ThreadLocal G4GlobalMagFieldMessengerfMagFieldMessenger = 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 class to define materials and geometry. The calorimeter is a box made of a given number of layers. A layer consists of an absorber plate and of a detection gap. The layer is replicated.

Four parameters define the geometry of the calorimeter :

  • the thickness of an absorber plate,
  • the thickness of a gap,
  • the number of layers,
  • the transverse size of the calorimeter (the input face is a square).

In ConstructSDandField() sensitive detectors of G4MultiFunctionalDetector type with primitive scorers are created and associated with the Absorber and Gap volumes. In addition a transverse uniform magnetic field is defined via G4GlobalMagFieldMessenger class.

Definition at line 56 of file B4dDetectorConstruction.hh.

Constructor & Destructor Documentation

◆ B4dDetectorConstruction()

B4dDetectorConstruction::B4dDetectorConstruction ( )

◆ ~B4dDetectorConstruction()

B4dDetectorConstruction::~B4dDetectorConstruction ( )
virtual

Definition at line 71 of file B4dDetectorConstruction.cc.

72 {
73 }

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * B4dDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 77 of file B4dDetectorConstruction.cc.

78 {
79  // Define materials
81 
82  // Define volumes
83  return DefineVolumes();
84 }
G4VPhysicalVolume * DefineVolumes()
Here is the call graph for this function:

◆ ConstructSDandField()

void B4dDetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 279 of file B4dDetectorConstruction.cc.

280 {
282  //
283  // Scorers
284  //
285 
286  // declare Absorber as a MultiFunctionalDetector scorer
287  //
288  G4MultiFunctionalDetector* absDetector
289  = new G4MultiFunctionalDetector("Absorber");
290 
291  G4VPrimitiveScorer* primitive;
292  primitive = new G4PSEnergyDeposit("Edep");
293  absDetector->RegisterPrimitive(primitive);
294 
295  primitive = new G4PSTrackLength("TrackLength");
296  G4SDChargedFilter* charged = new G4SDChargedFilter("chargedFilter");
297  primitive ->SetFilter(charged);
298  absDetector->RegisterPrimitive(primitive);
299 
300  SetSensitiveDetector("AbsoLV",absDetector);
301 
302  // declare Gap as a MultiFunctionalDetector scorer
303  //
304  G4MultiFunctionalDetector* gapDetector
305  = new G4MultiFunctionalDetector("Gap");
306 
307  primitive = new G4PSEnergyDeposit("Edep");
308  gapDetector->RegisterPrimitive(primitive);
309 
310  primitive = new G4PSTrackLength("TrackLength");
311  primitive ->SetFilter(charged);
312  gapDetector->RegisterPrimitive(primitive);
313 
314  SetSensitiveDetector("GapLV",gapDetector);
315 
316  //
317  // Magnetic field
318  //
319  // Create global magnetic field messenger.
320  // Uniform magnetic field is then created automatically if
321  // the field value is not zero.
322  G4ThreeVector fieldValue = G4ThreeVector();
325 
326  // Register the field messenger for deleting
328 }
void SetVerboseLevel(G4int verboseLevel)
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
CLHEP::Hep3Vector G4ThreeVector
void SetVerboseLevel(G4int vl)
Definition: G4SDManager.hh:90
void SetFilter(G4VSDFilter *f)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
static G4ThreadLocal G4GlobalMagFieldMessenger * fMagFieldMessenger
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
Here is the call graph for this function:

◆ DefineMaterials()

void B4dDetectorConstruction::DefineMaterials ( )
private

Definition at line 88 of file B4dDetectorConstruction.cc.

89 {
90  // Lead material defined using NIST Manager
91  G4NistManager* nistManager = G4NistManager::Instance();
92  nistManager->FindOrBuildMaterial("G4_Pb");
93 
94  // Liquid argon material
95  G4double a; // mass of a mole;
96  G4double z; // z=mean number of protons;
98  new G4Material("liquidArgon", z=18., a= 39.95*g/mole, density= 1.390*g/cm3);
99  // The argon by NIST Manager is a gas with a different density
100 
101  // Vacuum
102  new G4Material("Galactic", z=1., a=1.01*g/mole,density= universe_mean_density,
103  kStateGas, 2.73*kelvin, 3.e-18*pascal);
104 
105  // Print materials
107 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
int universe_mean_density
Definition: hepunit.py:307
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
static G4NistManager * Instance()
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
G4GLOB_DLL std::ostream G4cout
#define pascal
static const double cm3
Definition: G4SIunits.hh:120
static const double kelvin
Definition: G4SIunits.hh:278
static const double mole
Definition: G4SIunits.hh:283
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DefineVolumes()

G4VPhysicalVolume * B4dDetectorConstruction::DefineVolumes ( )
private

Definition at line 111 of file B4dDetectorConstruction.cc.

112 {
113  // Geometry parameters
114  G4int nofLayers = 10;
115  G4double absoThickness = 10.*mm;
116  G4double gapThickness = 5.*mm;
117  G4double calorSizeXY = 10.*cm;
118 
119  G4double layerThickness = absoThickness + gapThickness;
120  G4double calorThickness = nofLayers * layerThickness;
121  G4double worldSizeXY = 1.2 * calorSizeXY;
122  G4double worldSizeZ = 1.2 * calorThickness;
123 
124  // Get materials
125  G4Material* defaultMaterial = G4Material::GetMaterial("Galactic");
126  G4Material* absorberMaterial = G4Material::GetMaterial("G4_Pb");
127  G4Material* gapMaterial = G4Material::GetMaterial("liquidArgon");
128 
129  if ( ! defaultMaterial || ! absorberMaterial || ! gapMaterial ) {
131  msg << "Cannot retrieve materials already defined.";
132  G4Exception("B4DetectorConstruction::DefineVolumes()",
133  "MyCode0001", FatalException, msg);
134  }
135 
136  //
137  // World
138  //
139  G4VSolid* worldS
140  = new G4Box("World", // its name
141  worldSizeXY/2, worldSizeXY/2, worldSizeZ/2); // its size
142 
143  G4LogicalVolume* worldLV
144  = new G4LogicalVolume(
145  worldS, // its solid
146  defaultMaterial, // its material
147  "World"); // its name
148 
149  G4VPhysicalVolume* worldPV
150  = new G4PVPlacement(
151  0, // no rotation
152  G4ThreeVector(), // at (0,0,0)
153  worldLV, // its logical volume
154  "World", // its name
155  0, // its mother volume
156  false, // no boolean operation
157  0, // copy number
158  fCheckOverlaps); // checking overlaps
159 
160  //
161  // Calorimeter
162  //
163  G4VSolid* calorimeterS
164  = new G4Box("Calorimeter", // its name
165  calorSizeXY/2, calorSizeXY/2, calorThickness/2); // its size
166 
167  G4LogicalVolume* calorLV
168  = new G4LogicalVolume(
169  calorimeterS, // its solid
170  defaultMaterial, // its material
171  "Calorimeter"); // its name
172 
173  new G4PVPlacement(
174  0, // no rotation
175  G4ThreeVector(), // at (0,0,0)
176  calorLV, // its logical volume
177  "Calorimeter", // its name
178  worldLV, // its mother volume
179  false, // no boolean operation
180  0, // copy number
181  fCheckOverlaps); // checking overlaps
182 
183  //
184  // Layer
185  //
186  G4VSolid* layerS
187  = new G4Box("Layer", // its name
188  calorSizeXY/2, calorSizeXY/2, layerThickness/2); // its size
189 
190  G4LogicalVolume* layerLV
191  = new G4LogicalVolume(
192  layerS, // its solid
193  defaultMaterial, // its material
194  "Layer"); // its name
195 
196  new G4PVReplica(
197  "Layer", // its name
198  layerLV, // its logical volume
199  calorLV, // its mother
200  kZAxis, // axis of replication
201  nofLayers, // number of replica
202  layerThickness); // witdth of replica
203 
204  //
205  // Absorber
206  //
207  G4VSolid* absorberS
208  = new G4Box("Abso", // its name
209  calorSizeXY/2, calorSizeXY/2, absoThickness/2); // its size
210 
211  G4LogicalVolume* absorberLV
212  = new G4LogicalVolume(
213  absorberS, // its solid
214  absorberMaterial, // its material
215  "AbsoLV"); // its name
216 
217  new G4PVPlacement(
218  0, // no rotation
219  G4ThreeVector(0., 0., -gapThickness/2), // its position
220  absorberLV, // its logical volume
221  "Abso", // its name
222  layerLV, // its mother volume
223  false, // no boolean operation
224  0, // copy number
225  fCheckOverlaps); // checking overlaps
226 
227  //
228  // Gap
229  //
230  G4VSolid* gapS
231  = new G4Box("Gap", // its name
232  calorSizeXY/2, calorSizeXY/2, gapThickness/2); // its size
233 
234  G4LogicalVolume* gapLV
235  = new G4LogicalVolume(
236  gapS, // its solid
237  gapMaterial, // its material
238  "GapLV"); // its name
239 
240  new G4PVPlacement(
241  0, // no rotation
242  G4ThreeVector(0., 0., absoThickness/2), // its position
243  gapLV, // its logical volume
244  "Gap", // its name
245  layerLV, // its mother volume
246  false, // no boolean operation
247  0, // copy number
248  fCheckOverlaps); // checking overlaps
249 
250  //
251  // print parameters
252  //
253  G4cout
254  << G4endl
255  << "------------------------------------------------------------" << G4endl
256  << "---> The calorimeter is " << nofLayers << " layers of: [ "
257  << absoThickness/mm << "mm of " << absorberMaterial->GetName()
258  << " + "
259  << gapThickness/mm << "mm of " << gapMaterial->GetName() << " ] " << G4endl
260  << "------------------------------------------------------------" << G4endl;
261 
262  //
263  // Visualization attributes
264  //
266 
267  G4VisAttributes* simpleBoxVisAtt= new G4VisAttributes(G4Colour(1.0,1.0,1.0));
268  simpleBoxVisAtt->SetVisibility(true);
269  calorLV->SetVisAttributes(simpleBoxVisAtt);
270 
271  //
272  // Always return the physical World
273  //
274  return worldPV;
275 }
static const double cm
Definition: G4SIunits.hh:118
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Box.hh:64
void SetVisibility(G4bool)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4VisAttributes Invisible
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
static const double mm
Definition: G4SIunits.hh:114
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fCheckOverlaps

G4bool B4dDetectorConstruction::fCheckOverlaps
private

Definition at line 77 of file B4dDetectorConstruction.hh.

◆ fMagFieldMessenger

G4ThreadLocal G4GlobalMagFieldMessenger * B4dDetectorConstruction::fMagFieldMessenger = 0
staticprivate

Definition at line 74 of file B4dDetectorConstruction.hh.


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