Geant4  10.02.p03
ElectronBenchmarkDetector Class Reference

#include <ElectronBenchmarkDetector.hh>

Inheritance diagram for ElectronBenchmarkDetector:
Collaboration diagram for ElectronBenchmarkDetector:

Public Member Functions

 ElectronBenchmarkDetector ()
 
virtual ~ElectronBenchmarkDetector ()
 
virtual G4VPhysicalVolumeConstruct ()
 
void ConstructSDandField ()
 
void DefineMaterials ()
 
G4VPhysicalVolumeCreateWorld ()
 
void CreateExitWindow (G4LogicalVolume *logicWorld)
 
void CreatePrimaryFoil (G4LogicalVolume *logicWorld)
 
void CreateMonitor (G4LogicalVolume *logicWorld)
 
void CreateHeliumBag (G4LogicalVolume *logicWorld)
 
void CreateScorer (G4LogicalVolume *logicWorld)
 
G4VPhysicalVolumeCreateGeometry ()
 
void SetPrimFoilMaterial (G4String matname)
 
void SetPrimFoilThickness (G4double thicknessPrimFoil)
 
- 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 Attributes

G4double fPosWindow0
 
G4double fPosWindow1
 
G4MaterialfMaterialPrimFoil
 
G4double fHalfThicknessPrimFoil
 
G4double fPosPrimFoil
 
G4LogicalVolumefLogPrimFoil
 
G4TubsfSolidPrimFoil
 
G4double fPosMon0
 
G4double fPosMon1
 
G4double fPosBag0
 
G4double fPosBag1
 
G4double fPosHelium0
 
G4double fPosHelium1
 
G4double fThicknessRing
 
G4double fPosScorer
 
G4double fThicknessScorer
 
G4double fWidthScorerRing
 
G4LogicalVolumefScorerRingLog
 
G4double fRadOverall
 
G4double fRadRingInner
 
G4double fPosDelta
 
G4double fRadDelta
 
G4LogicalVolumefLogWorld
 
const G4Cache< G4MultiFunctionalDetector * > fSensitiveDetectorCache
 
ElectronBenchmarkDetectorMessengerfMessenger
 
G4VisAttributesfWorldVisAtt
 
G4VisAttributesfWindowVisAtt
 
G4VisAttributesfPrimFoilVisAtt
 
G4VisAttributesfMonVisAtt
 
G4VisAttributesfBagVisAtt
 
G4VisAttributesfHeliumVisAtt
 
G4VisAttributesfRingVisAtt
 
G4VisAttributesfScorerVisAtt
 

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

Definition at line 49 of file ElectronBenchmarkDetector.hh.

Constructor & Destructor Documentation

◆ ElectronBenchmarkDetector()

ElectronBenchmarkDetector::ElectronBenchmarkDetector ( )

Definition at line 59 of file ElectronBenchmarkDetector.cc.

62 fLogPrimFoil(0),
65 fLogWorld(0),
66 fMessenger(0),
67 fWorldVisAtt(0),
68 fWindowVisAtt(0),
70 fMonVisAtt(0),
71 fBagVisAtt(0),
72 fHeliumVisAtt(0),
73 fRingVisAtt(0),
75 {
76  // Exit Window
77  fPosWindow0 = 0.000000*cm;
78  fPosWindow1 = 0.004120*cm;
79 
80  // Scattering Foil
81  fPosPrimFoil = 2.650000*cm;
83 
84  // Monitor Chamber
85  fPosMon0 = 5.000000*cm;
86  fPosMon1 = 5.011270*cm;
87 
88  // Helium Bag
89  fPosBag0 = 6.497500*cm;
90  fPosHelium0 = 6.500000*cm;
91  fPosHelium1 = 116.500000*cm;
92  fPosBag1 = 116.502500*cm;
93  fThicknessRing = 1.4*cm;
94 
95  // Scoring Plane
96  fPosScorer = 118.200000*cm;
97  fThicknessScorer= 0.001*cm;
98  fWidthScorerRing= 0.1*cm;
99 
100  // Radii
101  fRadOverall = 23.3*cm;
102  fRadRingInner = 20.0*cm;
103 
104  // Extra space remaining in world volume around apparatus
105  fPosDelta = 1.*cm;
106  fRadDelta = 0.1*cm;
107 
109  DefineMaterials();
110 }
static const double cm
Definition: G4SIunits.hh:118
ElectronBenchmarkDetectorMessenger * fMessenger
Here is the call graph for this function:

◆ ~ElectronBenchmarkDetector()

ElectronBenchmarkDetector::~ElectronBenchmarkDetector ( )
virtual

Definition at line 114 of file ElectronBenchmarkDetector.cc.

115 {
116  delete fMessenger;
117 
118  delete fWorldVisAtt;
119  delete fWindowVisAtt;
120  delete fPrimFoilVisAtt;
121  delete fMonVisAtt;
122  delete fBagVisAtt;
123  delete fHeliumVisAtt;
124  delete fRingVisAtt;
125  delete fScorerVisAtt;
126 }
ElectronBenchmarkDetectorMessenger * fMessenger

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * ElectronBenchmarkDetector::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 130 of file ElectronBenchmarkDetector.cc.

131 {
132  return CreateGeometry();
133 }
Here is the call graph for this function:

◆ ConstructSDandField()

void ElectronBenchmarkDetector::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 390 of file ElectronBenchmarkDetector.cc.

391 {
393 
394  // G4Cache mechanism is necessary for multi-threaded operation
395  // as it allows us to store separate detector pointer per thread
396  G4MultiFunctionalDetector*& sensitiveDetector =
398 
399  if (!sensitiveDetector) {
400  sensitiveDetector = new G4MultiFunctionalDetector("MyDetector");
401 
402  G4VPrimitiveScorer* primitive;
403 
404  G4SDParticleFilter* electronFilter =
405  new G4SDParticleFilter("electronFilter", "e-");
406 
407  primitive = new G4PSCellFlux("cell flux");
408  sensitiveDetector->RegisterPrimitive(primitive);
409 
410  primitive = new G4PSCellFlux("e cell flux");
411  primitive->SetFilter(electronFilter);
412  sensitiveDetector->RegisterPrimitive(primitive);
413 
414  primitive = new G4PSPopulation("population");
415  sensitiveDetector->RegisterPrimitive(primitive);
416 
417  primitive = new G4PSPopulation("e population");
418  primitive->SetFilter(electronFilter);
419  sensitiveDetector->RegisterPrimitive(primitive);
420  }
421 
422  SetSensitiveDetector("scorerRingLog",sensitiveDetector);
423 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
const G4Cache< G4MultiFunctionalDetector * > fSensitiveDetectorCache
void SetVerboseLevel(G4int vl)
Definition: G4SDManager.hh:90
void SetFilter(G4VSDFilter *f)
value_type & Get() const
Definition: G4Cache.hh:282
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateExitWindow()

void ElectronBenchmarkDetector::CreateExitWindow ( G4LogicalVolume logicWorld)

Definition at line 235 of file ElectronBenchmarkDetector.cc.

235  {
236  G4double halfLengthWorld = fPosScorer/2.;
237  G4double halfThicknessWindow = fPosWindow1/2.;
238  G4VSolid* windowSolid = new G4Tubs("windowSolid", 0.*cm, fRadOverall,
239  halfThicknessWindow, 0.*deg, 360.*deg);
240  G4LogicalVolume* windowLog = new G4LogicalVolume(windowSolid,
241  G4Material::GetMaterial("TiAlloy"),
242  "windowLog");
243 
244  fWindowVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
245  windowLog->SetVisAttributes(fWindowVisAtt);
246 
247  new G4PVPlacement(0,
248  G4ThreeVector(0.,0.,
249  halfThicknessWindow - halfLengthWorld),
250  windowLog,"ExitWindow",worldLog,false,0);
251 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Tubs.hh:85
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateGeometry()

G4VPhysicalVolume * ElectronBenchmarkDetector::CreateGeometry ( )

Definition at line 190 of file ElectronBenchmarkDetector.cc.

190  {
191  // Clean old geometry, if any
196 
197  // Instantiate the world
198  G4VPhysicalVolume* physiworld = CreateWorld();
199  fLogWorld = physiworld->GetLogicalVolume();
200 
201  // Instantiate the geometry
206 
207  // Create the scorers
209 
210  return physiworld;
211 }
void CreateScorer(G4LogicalVolume *logicWorld)
static void Clean()
Definition: G4SolidStore.cc:79
void CreatePrimaryFoil(G4LogicalVolume *logicWorld)
static G4PhysicalVolumeStore * GetInstance()
void CreateHeliumBag(G4LogicalVolume *logicWorld)
void CreateMonitor(G4LogicalVolume *logicWorld)
static G4LogicalVolumeStore * GetInstance()
static G4SolidStore * GetInstance()
void CreateExitWindow(G4LogicalVolume *logicWorld)
static G4GeometryManager * GetInstance()
void OpenGeometry(G4VPhysicalVolume *vol=0)
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateHeliumBag()

void ElectronBenchmarkDetector::CreateHeliumBag ( G4LogicalVolume logicWorld)

Definition at line 297 of file ElectronBenchmarkDetector.cc.

297  {
298  G4double halfLengthWorld = fPosScorer/2.;
299 
300  // Construct cylinder of Mylar
301  G4double halfThicknessBag = (fPosBag1 - fPosBag0) /2.;
302  G4VSolid* bagSolid = new G4Tubs("bagSolid", 0.*cm, fRadOverall,
303  halfThicknessBag, 0.*deg, 360.*deg);
304  G4LogicalVolume* bagLog = new G4LogicalVolume(bagSolid,
305  G4Material::GetMaterial("G4_MYLAR"),
306  "bagLog");
307 
308  fBagVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
309  bagLog->SetVisAttributes(fBagVisAtt);
310 
311  new G4PVPlacement(0,
312  G4ThreeVector(0.,0.,
313  fPosBag0 + halfThicknessBag - halfLengthWorld),
314  bagLog,"HeliumBag",worldLog,false,0);
315 
316  // Insert cylinder of Helium into the Cylinder of Mylar
317  G4double halfThicknessHelium = (fPosHelium1 - fPosHelium0) /2.;
318  G4VSolid* heliumSolid = new G4Tubs("heliumSolid", 0.*cm, fRadOverall,
319  halfThicknessHelium, 0.*deg, 360.*deg);
320  G4LogicalVolume* heliumLog = new G4LogicalVolume(heliumSolid,
321  G4Material::GetMaterial("G4_He"),
322  "heliumLog");
323 
324  fHeliumVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
325  heliumLog->SetVisAttributes(fHeliumVisAtt);
326 
327  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.),
328  heliumLog,"Helium",bagLog,false,0);
329 
330  // Insert two rings of Aluminum into the Cylinder of Helium
331  G4double halfThicknessRing = fThicknessRing /2.;
332  G4VSolid* ringSolid = new G4Tubs("ringSolid", fRadRingInner, fRadOverall,
333  halfThicknessRing, 0.*deg, 360.*deg);
334  G4LogicalVolume* ring0Log = new G4LogicalVolume(ringSolid,
335  G4Material::GetMaterial("G4_Al"),
336  "ring0Log");
337  G4LogicalVolume* ring1Log = new G4LogicalVolume(ringSolid,
338  G4Material::GetMaterial("G4_Al"),
339  "ring1Log");
340 
341  fRingVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
342  ring0Log->SetVisAttributes(fRingVisAtt);
343  ring1Log->SetVisAttributes(fRingVisAtt);
344 
345  new G4PVPlacement(0,
346  G4ThreeVector(0.,0.,
347  -halfThicknessHelium + halfThicknessRing),
348  ring0Log,"Ring0",heliumLog,false,0);
349 
350  new G4PVPlacement(0,
351  G4ThreeVector(0.,0.,
352  halfThicknessHelium - halfThicknessRing),
353  ring1Log,"Ring1",heliumLog,false,0);
354 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Tubs.hh:85
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateMonitor()

void ElectronBenchmarkDetector::CreateMonitor ( G4LogicalVolume logicWorld)

Definition at line 277 of file ElectronBenchmarkDetector.cc.

277  {
278  G4double halfLengthWorld = fPosScorer/2.;
279  G4double halfThicknessMon = (fPosMon1 - fPosMon0) /2.;
280  G4VSolid* monSolid = new G4Tubs("monSolid", 0.*cm, fRadOverall,
281  halfThicknessMon, 0.*deg, 360.*deg);
282  G4LogicalVolume* monLog = new G4LogicalVolume(monSolid,
283  G4Material::GetMaterial("G4_MYLAR"),
284  "monLog");
285 
286  fMonVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
287  monLog->SetVisAttributes(fMonVisAtt);
288 
289  new G4PVPlacement(0,
290  G4ThreeVector(0.,0.,
291  fPosMon0 + halfThicknessMon - halfLengthWorld),
292  monLog,"MonitorChamber",worldLog,false,0);
293 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Tubs.hh:85
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreatePrimaryFoil()

void ElectronBenchmarkDetector::CreatePrimaryFoil ( G4LogicalVolume logicWorld)

Definition at line 255 of file ElectronBenchmarkDetector.cc.

255  {
256  G4double halfLengthWorld = fPosScorer/2.;
257 
258  // For some energies, we have no Primary Foil.
259  if (fHalfThicknessPrimFoil==0.) return;
260 
261  fSolidPrimFoil = new G4Tubs("PrimFoilSolid", 0.*cm, fRadOverall,
262  fHalfThicknessPrimFoil, 0.*deg, 360.*deg);
264  fMaterialPrimFoil, "PrimFoilLog");
265 
266  fPrimFoilVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
268 
269  new G4PVPlacement(0,
270  G4ThreeVector(0.,0.,
271  fPosPrimFoil + fHalfThicknessPrimFoil - halfLengthWorld),
272  fLogPrimFoil,"ScatteringFoil",worldLog,false,0);
273 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Tubs.hh:85
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateScorer()

void ElectronBenchmarkDetector::CreateScorer ( G4LogicalVolume logicWorld)

Definition at line 358 of file ElectronBenchmarkDetector.cc.

358  {
359  G4double halfLengthWorld = fPosScorer/2.;
360  G4double halfThicknessScorer = fThicknessScorer /2.;
361 
362  G4VSolid* scorerSolid = new G4Tubs("scorerSolid", 0.*cm, fRadOverall,
363  halfThicknessScorer, 0.*deg, 360.*deg);
364  G4LogicalVolume* scorerLog = new G4LogicalVolume(scorerSolid,
365  G4Material::GetMaterial("G4_AIR"),
366  "scorerLog");
367 
368  fScorerVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
369  scorerLog->SetVisAttributes(fScorerVisAtt);
370  new G4PVPlacement(0,
371  G4ThreeVector(0.,0.,
372  halfLengthWorld - halfThicknessScorer),
373  scorerLog,"Scorer",worldLog,false,0);
374 
375  G4VSolid* scorerRingSolid = new G4Tubs("scorerRingSolid", 0.*cm,
376  fRadOverall,
377  halfThicknessScorer, 0.*deg, 360.*deg);
378  fScorerRingLog = new G4LogicalVolume(scorerRingSolid,
379  G4Material::GetMaterial("G4_AIR"), "scorerRingLog");
380  new G4PVReplica("ScorerRing",fScorerRingLog,scorerLog,kRho,
382 
384 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Tubs.hh:85
int G4int
Definition: G4Types.hh:78
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateWorld()

G4VPhysicalVolume * ElectronBenchmarkDetector::CreateWorld ( )

Definition at line 215 of file ElectronBenchmarkDetector.cc.

215  {
216  G4double halfLengthWorld = fPosScorer/2. + fPosDelta;
217  G4double radWorld = fRadOverall + fRadDelta;
218  G4VSolid* worldSolid = new G4Tubs("WorldSolid", 0.*cm, radWorld,
219  halfLengthWorld, 0.*deg, 360.*deg);
220  G4LogicalVolume* worldLog = new G4LogicalVolume(worldSolid,
221  G4Material::GetMaterial("G4_AIR"), "WorldLog");
222 
223  fWorldVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
224  worldLog->SetVisAttributes(fWorldVisAtt);
225 
226  G4VPhysicalVolume* worldPhys =
227  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.),
228  worldLog,"World", 0, false, 0);
229 
230  return worldPhys;
231 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
Definition: G4Tubs.hh:85
static const double deg
Definition: G4SIunits.hh:151
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DefineMaterials()

void ElectronBenchmarkDetector::DefineMaterials ( void  )

Definition at line 137 of file ElectronBenchmarkDetector.cc.

137  {
138  // Use NIST database for elements and materials whereever possible.
140  man->SetVerbose(1);
141 
142  // Take all elements and materials from NIST
143  man->FindOrBuildMaterial("G4_He");
144  man->FindOrBuildMaterial("G4_Be");
145  man->FindOrBuildMaterial("G4_Al");
146  man->FindOrBuildMaterial("G4_Ti");
147  man->FindOrBuildMaterial("G4_Ta");
148  man->FindOrBuildMaterial("G4_AIR");
149  man->FindOrBuildMaterial("G4_MYLAR");
150 
151  G4Element* C = man->FindOrBuildElement("C");
152  G4Element* Cu = man->FindOrBuildElement("Cu");
153  G4Element* Au = man->FindOrBuildElement("Au");
154  G4Element* Ti = man->FindOrBuildElement("Ti");
155  G4Element* Al = man->FindOrBuildElement("Al");
156  G4Element* V = man->FindOrBuildElement("V");
157 
158  // Define materials not in NIST.
159  // While the NIST database does contain default materials for C, Cu and Au,
160  // those defaults have different densities than the ones used in the
161  // benchmark specification.
163  G4int ncomponents;
164  G4double fractionmass;
165 
166  G4Material* G4_C = new G4Material("G4_C", density= 2.18*g/cm3,
167  ncomponents=1);
168  G4_C->AddElement(C, fractionmass=1.00);
169 
170  G4Material* G4_Cu = new G4Material("G4_Cu", density= 8.92*g/cm3,
171  ncomponents=1);
172  G4_Cu->AddElement(Cu, fractionmass=1.00);
173 
174  G4Material* G4_Au = new G4Material("G4_Au", density= 19.30*g/cm3,
175  ncomponents=1);
176  G4_Au->AddElement(Au, fractionmass=1.00);
177 
178  G4Material* TiAlloy = new G4Material("TiAlloy", density= 4.42*g/cm3,
179  ncomponents=3);
180  TiAlloy->AddElement(Ti, fractionmass=0.90);
181  TiAlloy->AddElement(Al, fractionmass=0.06);
182  TiAlloy->AddElement(V, fractionmass=0.04);
183 
184  // Print materials table
186 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
double C(double temp)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
G4GLOB_DLL std::ostream G4cout
void SetVerbose(G4int)
static const double cm3
Definition: G4SIunits.hh:120
#define G4endl
Definition: G4ios.hh:61
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:364
double G4double
Definition: G4Types.hh:76
G4Element * FindOrBuildElement(G4int Z, G4bool isotopes=true)
G4Material * Al
Definition: TRTMaterials.hh:74
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPrimFoilMaterial()

void ElectronBenchmarkDetector::SetPrimFoilMaterial ( G4String  matname)

Definition at line 427 of file ElectronBenchmarkDetector.cc.

427  {
429  if (fLogPrimFoil) {
431  }
434 }
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
void CreatePrimaryFoil(G4LogicalVolume *logicWorld)
void PhysicsHasBeenModified()
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
void SetMaterial(G4Material *pMaterial)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPrimFoilThickness()

void ElectronBenchmarkDetector::SetPrimFoilThickness ( G4double  thicknessPrimFoil)

Definition at line 438 of file ElectronBenchmarkDetector.cc.

439 {
440  fHalfThicknessPrimFoil = thicknessPrimFoil / 2.;
441  if (fSolidPrimFoil) {
443  }
446 }
void GeometryHasBeenModified(G4bool prop=true)
void CreatePrimaryFoil(G4LogicalVolume *logicWorld)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
void SetZHalfLength(G4double newDz)
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fBagVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fBagVisAtt
private

Definition at line 130 of file ElectronBenchmarkDetector.hh.

◆ fHalfThicknessPrimFoil

G4double ElectronBenchmarkDetector::fHalfThicknessPrimFoil
private

Definition at line 84 of file ElectronBenchmarkDetector.hh.

◆ fHeliumVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fHeliumVisAtt
private

Definition at line 131 of file ElectronBenchmarkDetector.hh.

◆ fLogPrimFoil

G4LogicalVolume* ElectronBenchmarkDetector::fLogPrimFoil
private

Definition at line 86 of file ElectronBenchmarkDetector.hh.

◆ fLogWorld

G4LogicalVolume* ElectronBenchmarkDetector::fLogWorld
private

Definition at line 115 of file ElectronBenchmarkDetector.hh.

◆ fMaterialPrimFoil

G4Material* ElectronBenchmarkDetector::fMaterialPrimFoil
private

Definition at line 83 of file ElectronBenchmarkDetector.hh.

◆ fMessenger

ElectronBenchmarkDetectorMessenger* ElectronBenchmarkDetector::fMessenger
private

Definition at line 123 of file ElectronBenchmarkDetector.hh.

◆ fMonVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fMonVisAtt
private

Definition at line 129 of file ElectronBenchmarkDetector.hh.

◆ fPosBag0

G4double ElectronBenchmarkDetector::fPosBag0
private

Definition at line 94 of file ElectronBenchmarkDetector.hh.

◆ fPosBag1

G4double ElectronBenchmarkDetector::fPosBag1
private

Definition at line 95 of file ElectronBenchmarkDetector.hh.

◆ fPosDelta

G4double ElectronBenchmarkDetector::fPosDelta
private

Definition at line 111 of file ElectronBenchmarkDetector.hh.

◆ fPosHelium0

G4double ElectronBenchmarkDetector::fPosHelium0
private

Definition at line 96 of file ElectronBenchmarkDetector.hh.

◆ fPosHelium1

G4double ElectronBenchmarkDetector::fPosHelium1
private

Definition at line 97 of file ElectronBenchmarkDetector.hh.

◆ fPosMon0

G4double ElectronBenchmarkDetector::fPosMon0
private

Definition at line 90 of file ElectronBenchmarkDetector.hh.

◆ fPosMon1

G4double ElectronBenchmarkDetector::fPosMon1
private

Definition at line 91 of file ElectronBenchmarkDetector.hh.

◆ fPosPrimFoil

G4double ElectronBenchmarkDetector::fPosPrimFoil
private

Definition at line 85 of file ElectronBenchmarkDetector.hh.

◆ fPosScorer

G4double ElectronBenchmarkDetector::fPosScorer
private

Definition at line 101 of file ElectronBenchmarkDetector.hh.

◆ fPosWindow0

G4double ElectronBenchmarkDetector::fPosWindow0
private

Definition at line 79 of file ElectronBenchmarkDetector.hh.

◆ fPosWindow1

G4double ElectronBenchmarkDetector::fPosWindow1
private

Definition at line 80 of file ElectronBenchmarkDetector.hh.

◆ fPrimFoilVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fPrimFoilVisAtt
private

Definition at line 128 of file ElectronBenchmarkDetector.hh.

◆ fRadDelta

G4double ElectronBenchmarkDetector::fRadDelta
private

Definition at line 112 of file ElectronBenchmarkDetector.hh.

◆ fRadOverall

G4double ElectronBenchmarkDetector::fRadOverall
private

Definition at line 107 of file ElectronBenchmarkDetector.hh.

◆ fRadRingInner

G4double ElectronBenchmarkDetector::fRadRingInner
private

Definition at line 108 of file ElectronBenchmarkDetector.hh.

◆ fRingVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fRingVisAtt
private

Definition at line 132 of file ElectronBenchmarkDetector.hh.

◆ fScorerRingLog

G4LogicalVolume* ElectronBenchmarkDetector::fScorerRingLog
private

Definition at line 104 of file ElectronBenchmarkDetector.hh.

◆ fScorerVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fScorerVisAtt
private

Definition at line 133 of file ElectronBenchmarkDetector.hh.

◆ fSensitiveDetectorCache

const G4Cache<G4MultiFunctionalDetector*> ElectronBenchmarkDetector::fSensitiveDetectorCache
private

Definition at line 120 of file ElectronBenchmarkDetector.hh.

◆ fSolidPrimFoil

G4Tubs* ElectronBenchmarkDetector::fSolidPrimFoil
private

Definition at line 87 of file ElectronBenchmarkDetector.hh.

◆ fThicknessRing

G4double ElectronBenchmarkDetector::fThicknessRing
private

Definition at line 98 of file ElectronBenchmarkDetector.hh.

◆ fThicknessScorer

G4double ElectronBenchmarkDetector::fThicknessScorer
private

Definition at line 102 of file ElectronBenchmarkDetector.hh.

◆ fWidthScorerRing

G4double ElectronBenchmarkDetector::fWidthScorerRing
private

Definition at line 103 of file ElectronBenchmarkDetector.hh.

◆ fWindowVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fWindowVisAtt
private

Definition at line 127 of file ElectronBenchmarkDetector.hh.

◆ fWorldVisAtt

G4VisAttributes* ElectronBenchmarkDetector::fWorldVisAtt
private

Definition at line 126 of file ElectronBenchmarkDetector.hh.


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