Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
ElectronBenchmarkDetector Class Reference

#include <ElectronBenchmarkDetector.hh>

Inheritance diagram for ElectronBenchmarkDetector:
G4VUserDetectorConstruction

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 UpdateGeometry ()
 
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
 

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 48 of file ElectronBenchmarkDetector.hh.

Constructor & Destructor Documentation

ElectronBenchmarkDetector::ElectronBenchmarkDetector ( )

Definition at line 59 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, and DefineMaterials().

61 fMaterialPrimFoil(0),
62 fScorerRingLog(0),
63 fLogWorld(0),
64 fMessenger(0),
65 fWorldVisAtt(0),
66 fWindowVisAtt(0),
67 fPrimFoilVisAtt(0),
68 fMonVisAtt(0),
69 fBagVisAtt(0),
70 fHeliumVisAtt(0),
71 fRingVisAtt(0),
72 fScorerVisAtt(0)
73 {
74  // Exit Window
75  fPosWindow0 = 0.000000*cm;
76  fPosWindow1 = 0.004120*cm;
77 
78  // Scattering Foil
79  fPosPrimFoil = 2.650000*cm;
80 
81  // Monitor Chamber
82  fPosMon0 = 5.000000*cm;
83  fPosMon1 = 5.011270*cm;
84 
85  // Helium Bag
86  fPosBag0 = 6.497500*cm;
87  fPosHelium0 = 6.500000*cm;
88  fPosHelium1 = 116.500000*cm;
89  fPosBag1 = 116.502500*cm;
90  fThicknessRing = 1.4*cm;
91 
92  // Scoring Plane
93  fPosScorer = 118.200000*cm;
94  fThicknessScorer= 0.001*cm;
95  fWidthScorerRing= 0.1*cm;
96 
97  // Radii
98  fRadOverall = 23.3*cm;
99  fRadRingInner = 20.0*cm;
100 
101  // Extra space remaining in world volume around apparatus
102  fPosDelta = 1.*cm;
103  fRadDelta = 0.1*cm;
104 
105  fMessenger = new ElectronBenchmarkDetectorMessenger(this);
106  DefineMaterials();
107 }
ElectronBenchmarkDetector::~ElectronBenchmarkDetector ( )
virtual

Definition at line 111 of file ElectronBenchmarkDetector.cc.

112 {
113  delete fMessenger;
114 
115  delete fWorldVisAtt;
116  delete fWindowVisAtt;
117  delete fPrimFoilVisAtt;
118  delete fMonVisAtt;
119  delete fBagVisAtt;
120  delete fHeliumVisAtt;
121  delete fRingVisAtt;
122  delete fScorerVisAtt;
123 }

Member Function Documentation

G4VPhysicalVolume * ElectronBenchmarkDetector::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 127 of file ElectronBenchmarkDetector.cc.

References CreateGeometry().

128 {
129  return CreateGeometry();
130 }
void ElectronBenchmarkDetector::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 393 of file ElectronBenchmarkDetector.cc.

References G4Cache< VALTYPE >::Get(), G4SDManager::GetSDMpointer(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VPrimitiveScorer::SetFilter(), G4VUserDetectorConstruction::SetSensitiveDetector(), and G4SDManager::SetVerboseLevel().

Referenced by CreateScorer().

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

Definition at line 238 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, python.hepunit::deg, G4Material::GetMaterial(), and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

238  {
239  G4double halfLengthWorld = fPosScorer/2.;
240  G4double halfThicknessWindow = fPosWindow1/2.;
241  G4VSolid* windowSolid = new G4Tubs("windowSolid", 0.*cm, fRadOverall,
242  halfThicknessWindow, 0.*deg, 360.*deg);
243  G4LogicalVolume* windowLog = new G4LogicalVolume(windowSolid,
244  G4Material::GetMaterial("TiAlloy"),
245  "windowLog");
246 
247  fWindowVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
248  windowLog->SetVisAttributes(fWindowVisAtt);
249 
250  new G4PVPlacement(0,
251  G4ThreeVector(0.,0.,
252  halfThicknessWindow - halfLengthWorld),
253  windowLog,"ExitWindow",worldLog,false,0);
254 }
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
Definition: G4Tubs.hh:84
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
G4VPhysicalVolume * ElectronBenchmarkDetector::CreateGeometry ( )

Definition at line 187 of file ElectronBenchmarkDetector.cc.

References G4PhysicalVolumeStore::Clean(), G4SolidStore::Clean(), G4LogicalVolumeStore::Clean(), CreateExitWindow(), CreateHeliumBag(), CreateMonitor(), CreatePrimaryFoil(), CreateScorer(), CreateWorld(), G4SolidStore::GetInstance(), G4PhysicalVolumeStore::GetInstance(), G4LogicalVolumeStore::GetInstance(), G4GeometryManager::GetInstance(), G4VPhysicalVolume::GetLogicalVolume(), and G4GeometryManager::OpenGeometry().

Referenced by Construct().

187  {
188  // Clean old geometry, if any
193 
194  // Instantiate the world
195  G4VPhysicalVolume* physiworld = CreateWorld();
196  fLogWorld = physiworld->GetLogicalVolume();
197 
198  // Instantiate the geometry
199  CreateExitWindow(fLogWorld);
200  CreatePrimaryFoil(fLogWorld);
201  CreateMonitor(fLogWorld);
202  CreateHeliumBag(fLogWorld);
203 
204  // Create the scorers
205  CreateScorer(fLogWorld);
206 
207  return physiworld;
208 }
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()
G4LogicalVolume * GetLogicalVolume() const
void OpenGeometry(G4VPhysicalVolume *vol=0)
void ElectronBenchmarkDetector::CreateHeliumBag ( G4LogicalVolume logicWorld)

Definition at line 300 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, python.hepunit::deg, G4Material::GetMaterial(), and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

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

Definition at line 280 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, python.hepunit::deg, G4Material::GetMaterial(), and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

280  {
281  G4double halfLengthWorld = fPosScorer/2.;
282  G4double halfThicknessMon = (fPosMon1 - fPosMon0) /2.;
283  G4VSolid* monSolid = new G4Tubs("monSolid", 0.*cm, fRadOverall,
284  halfThicknessMon, 0.*deg, 360.*deg);
285  G4LogicalVolume* monLog = new G4LogicalVolume(monSolid,
286  G4Material::GetMaterial("G4_MYLAR"),
287  "monLog");
288 
289  fMonVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
290  monLog->SetVisAttributes(fMonVisAtt);
291 
292  new G4PVPlacement(0,
293  G4ThreeVector(0.,0.,
294  fPosMon0 + halfThicknessMon - halfLengthWorld),
295  monLog,"MonitorChamber",worldLog,false,0);
296 }
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
Definition: G4Tubs.hh:84
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
void ElectronBenchmarkDetector::CreatePrimaryFoil ( G4LogicalVolume logicWorld)

Definition at line 258 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, python.hepunit::deg, and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

258  {
259  G4double halfLengthWorld = fPosScorer/2.;
260 
261  // For some energies, we have no Primary Foil.
262  if (fHalfThicknessPrimFoil==0.) return;
263 
264  G4VSolid* primFoilSolid = new G4Tubs("PrimFoilSolid", 0.*cm, fRadOverall,
265  fHalfThicknessPrimFoil, 0.*deg, 360.*deg);
266  G4LogicalVolume* primFoilLog = new G4LogicalVolume(primFoilSolid,
267  fMaterialPrimFoil, "PrimFoilLog");
268 
269  fPrimFoilVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
270  primFoilLog->SetVisAttributes(fPrimFoilVisAtt);
271 
272  new G4PVPlacement(0,
273  G4ThreeVector(0.,0.,
274  fPosPrimFoil + fHalfThicknessPrimFoil - halfLengthWorld),
275  primFoilLog,"ScatteringFoil",worldLog,false,0);
276 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Tubs.hh:84
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
void ElectronBenchmarkDetector::CreateScorer ( G4LogicalVolume logicWorld)

Definition at line 361 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, ConstructSDandField(), python.hepunit::deg, G4Material::GetMaterial(), kRho, and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

361  {
362  G4double halfLengthWorld = fPosScorer/2.;
363  G4double halfThicknessScorer = fThicknessScorer /2.;
364 
365  G4VSolid* scorerSolid = new G4Tubs("scorerSolid", 0.*cm, fRadOverall,
366  halfThicknessScorer, 0.*deg, 360.*deg);
367  G4LogicalVolume* scorerLog = new G4LogicalVolume(scorerSolid,
368  G4Material::GetMaterial("G4_AIR"),
369  "scorerLog");
370 
371  fScorerVisAtt = new G4VisAttributes(G4Colour(0.5,1.0,0.5));
372  scorerLog->SetVisAttributes(fScorerVisAtt);
373  new G4PVPlacement(0,
374  G4ThreeVector(0.,0.,
375  halfLengthWorld - halfThicknessScorer),
376  scorerLog,"Scorer",worldLog,false,0);
377 
378  G4VSolid* scorerRingSolid = new G4Tubs("scorerRingSolid", 0.*cm,
379  fRadOverall,
380  halfThicknessScorer, 0.*deg, 360.*deg);
381  fScorerRingLog = new G4LogicalVolume(scorerRingSolid,
382  G4Material::GetMaterial("G4_AIR"), "scorerRingLog");
383  new G4PVReplica("ScorerRing",fScorerRingLog,scorerLog,kRho,
384  G4int(fRadOverall/fWidthScorerRing),fWidthScorerRing);
385 
387 }
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
Definition: G4Tubs.hh:84
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
void SetVisAttributes(const G4VisAttributes *pVA)
G4VPhysicalVolume * ElectronBenchmarkDetector::CreateWorld ( )

Definition at line 218 of file ElectronBenchmarkDetector.cc.

References python.hepunit::cm, python.hepunit::deg, G4Material::GetMaterial(), and G4LogicalVolume::SetVisAttributes().

Referenced by CreateGeometry().

218  {
219  G4double halfLengthWorld = fPosScorer/2. + fPosDelta;
220  G4double radWorld = fRadOverall + fRadDelta;
221  G4VSolid* worldSolid = new G4Tubs("WorldSolid", 0.*cm, radWorld,
222  halfLengthWorld, 0.*deg, 360.*deg);
223  G4LogicalVolume* worldLog = new G4LogicalVolume(worldSolid,
224  G4Material::GetMaterial("G4_AIR"), "WorldLog");
225 
226  fWorldVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
227  worldLog->SetVisAttributes(fWorldVisAtt);
228 
229  G4VPhysicalVolume* worldPhys =
230  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.),
231  worldLog,"World", 0, false, 0);
232 
233  return worldPhys;
234 }
CLHEP::Hep3Vector G4ThreeVector
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
Definition: G4Tubs.hh:84
double G4double
Definition: G4Types.hh:76
void SetVisAttributes(const G4VisAttributes *pVA)
void ElectronBenchmarkDetector::DefineMaterials ( )

Definition at line 134 of file ElectronBenchmarkDetector.cc.

References G4Material::AddElement(), Al, python.hepunit::cm3, density, G4NistManager::FindOrBuildElement(), G4NistManager::FindOrBuildMaterial(), g(), G4cout, G4endl, G4Material::GetMaterialTable(), G4NistManager::Instance(), and G4NistManager::SetVerbose().

Referenced by ElectronBenchmarkDetector().

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

Definition at line 430 of file ElectronBenchmarkDetector.cc.

References G4Material::GetMaterial(), and UpdateGeometry().

Referenced by ElectronBenchmarkDetectorMessenger::SetNewValue().

430  {
431  fMaterialPrimFoil = G4Material::GetMaterial(matname);
432  UpdateGeometry();
433 }
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
void ElectronBenchmarkDetector::SetPrimFoilThickness ( G4double  thicknessPrimFoil)

Definition at line 437 of file ElectronBenchmarkDetector.cc.

References UpdateGeometry().

Referenced by ElectronBenchmarkDetectorMessenger::SetNewValue().

438 {
439  fHalfThicknessPrimFoil = thicknessPrimFoil / 2.;
440  UpdateGeometry();
441 }
void ElectronBenchmarkDetector::UpdateGeometry ( )

Definition at line 212 of file ElectronBenchmarkDetector.cc.

References G4RunManager::GetRunManager(), and G4RunManager::ReinitializeGeometry().

Referenced by SetPrimFoilMaterial(), and SetPrimFoilThickness().

212  {
214 }
void ReinitializeGeometry(G4bool destroyFirst=false, G4bool prop=true)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74

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