Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B01DetectorConstruction Class Reference

#include <B01DetectorConstruction.hh>

Inheritance diagram for B01DetectorConstruction:
Collaboration diagram for B01DetectorConstruction:

Public Member Functions

 B01DetectorConstruction ()
 
 ~B01DetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
G4VIStoreCreateImportanceStore ()
 
G4VWeightWindowStoreCreateWeightWindowStore ()
 
G4String GetCellName (G4int i)
 
G4VPhysicalVolumeGetWorldVolume ()
 
void SetSensitive ()
 
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
 

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 44 of file B01DetectorConstruction.hh.

Constructor & Destructor Documentation

B01DetectorConstruction::B01DetectorConstruction ( )

Definition at line 68 of file B01DetectorConstruction.cc.

68  :
70  fLogicalVolumeVector(),fPhysicalVolumeVector()
71 {;}
B01DetectorConstruction::~B01DetectorConstruction ( )

Definition at line 75 of file B01DetectorConstruction.cc.

76 {
77  fLogicalVolumeVector.clear();
78  fPhysicalVolumeVector.clear();
79 }

Member Function Documentation

G4VPhysicalVolume * B01DetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 83 of file B01DetectorConstruction.cc.

84 {
85  G4double pos_x;
86  G4double pos_y;
87  G4double pos_z;
88 
89  G4double density, pressure, temperature;
90  G4double A;
91  G4int Z;
92 
93  G4String name, symbol;
94  G4double z;
95  G4double fractionmass;
96 
97  A = 1.01*g/mole;
98  G4Element* elH = new G4Element(name="Hydrogen",symbol="H" , Z= 1, A);
99 
100  A = 12.01*g/mole;
101  G4Element* elC = new G4Element(name="Carbon" ,symbol="C" , Z = 6, A);
102 
103  A = 16.00*g/mole;
104  G4Element* elO = new G4Element(name="Oxygen" ,symbol="O" , Z= 8, A);
105 
106  A = 22.99*g/mole;
107  G4Element* elNa = new G4Element(name="Natrium" ,symbol="Na" , Z=11 , A);
108 
109  A = 200.59*g/mole;
110  G4Element* elHg = new G4Element(name="Hg" ,symbol="Hg" , Z=80, A);
111 
112  A = 26.98*g/mole;
113  G4Element* elAl = new G4Element(name="Aluminium" ,symbol="Al" , Z=13, A);
114 
115  A = 28.09*g/mole;
116  G4Element* elSi = new G4Element(name="Silicon", symbol="Si", Z=14, A);
117 
118  A = 39.1*g/mole;
119  G4Element* elK = new G4Element(name="K" ,symbol="K" , Z=19 , A);
120 
121  A = 69.72*g/mole;
122  G4Element* elCa = new G4Element(name="Calzium" ,symbol="Ca" , Z=31 , A);
123 
124  A = 55.85*g/mole;
125  G4Element* elFe = new G4Element(name="Iron" ,symbol="Fe", Z=26, A);
126 
127  density = universe_mean_density; //from PhysicalConstants.h
128  pressure = 3.e-18*pascal;
129  temperature = 2.73*kelvin;
130  G4Material *Galactic =
131  new G4Material(name="Galactic", z=1., A=1.01*g/mole, density,
132  kStateGas,temperature,pressure);
133 
134  density = 2.03*g/cm3;
135  G4Material* Concrete = new G4Material("Concrete", density, 10);
136  Concrete->AddElement(elH , fractionmass= 0.01);
137  Concrete->AddElement(elO , fractionmass= 0.529);
138  Concrete->AddElement(elNa , fractionmass= 0.016);
139  Concrete->AddElement(elHg , fractionmass= 0.002);
140  Concrete->AddElement(elAl , fractionmass= 0.034);
141  Concrete->AddElement(elSi , fractionmass= 0.337);
142  Concrete->AddElement(elK , fractionmass= 0.013);
143  Concrete->AddElement(elCa , fractionmass= 0.044);
144  Concrete->AddElement(elFe , fractionmass= 0.014);
145  Concrete->AddElement(elC , fractionmass= 0.001);
146 
148  // world cylinder volume
150 
151  // world solid
152 
153  G4double innerRadiusCylinder = 0*cm;
154  G4double outerRadiusCylinder = 100*cm;
155  G4double heightCylinder = 100*cm;
156  G4double startAngleCylinder = 0*deg;
157  G4double spanningAngleCylinder = 360*deg;
158 
159  G4Tubs *worldCylinder = new G4Tubs("worldCylinder",
160  innerRadiusCylinder,
161  outerRadiusCylinder,
162  heightCylinder,
163  startAngleCylinder,
164  spanningAngleCylinder);
165 
166  // logical world
167 
168  G4LogicalVolume *worldCylinder_log =
169  new G4LogicalVolume(worldCylinder, Galactic, "worldCylinder_log");
170  fLogicalVolumeVector.push_back(worldCylinder_log);
171 
172  name = "shieldWorld";
173  fWorldVolume = new
174  G4PVPlacement(0, G4ThreeVector(0,0,0), worldCylinder_log,
175  name, 0, false, 0);
176 
177  fPhysicalVolumeVector.push_back(fWorldVolume);
178 
179  // creating 18 slabs of 10 cm thick concrete
180 
181  G4double innerRadiusShield = 0*cm;
182  G4double outerRadiusShield = 100*cm;
183  G4double heightShield = 5*cm;
184  G4double startAngleShield = 0*deg;
185  G4double spanningAngleShield = 360*deg;
186 
187  G4Tubs *aShield = new G4Tubs("aShield",
188  innerRadiusShield,
189  outerRadiusShield,
190  heightShield,
191  startAngleShield,
192  spanningAngleShield);
193 
194  // logical shield
195 
196  G4LogicalVolume *aShield_log =
197  new G4LogicalVolume(aShield, Concrete, "aShield_log");
198  fLogicalVolumeVector.push_back(aShield_log);
199 
200  G4VisAttributes* pShieldVis = new G4VisAttributes(G4Colour(0.0,0.0,1.0));
201  pShieldVis->SetForceSolid(true);
202  aShield_log->SetVisAttributes(pShieldVis);
203 
204  // physical shields
205 
206  G4int i;
207  G4double startz = -85*cm;
208  for (i=1; i<=18; i++)
209  {
210  name = GetCellName(i);
211  pos_x = 0*cm;
212  pos_y = 0*cm;
213  pos_z = startz + (i-1) * (2*heightShield);
214  G4VPhysicalVolume *pvol =
215  new G4PVPlacement(0,
216  G4ThreeVector(pos_x, pos_y, pos_z),
217  aShield_log,
218  name,
219  worldCylinder_log,
220  false,
221  i);
222  fPhysicalVolumeVector.push_back(pvol);
223  }
224 
225  // filling the rest of the world volume behind the concrete with
226  // another slab which should get the same importance value
227  // or lower weight bound as the last slab
228  //
229  innerRadiusShield = 0*cm;
230  outerRadiusShield = 100*cm;
231  heightShield = 5*cm;
232  startAngleShield = 0*deg;
233  spanningAngleShield = 360*deg;
234 
235  G4Tubs *aRest = new G4Tubs("Rest",
236  innerRadiusShield,
237  outerRadiusShield,
238  heightShield,
239  startAngleShield,
240  spanningAngleShield);
241 
242  G4LogicalVolume *aRest_log =
243  new G4LogicalVolume(aRest, Galactic, "aRest_log");
244  fLogicalVolumeVector.push_back(aRest_log);
245  name = "rest";
246 
247  pos_x = 0*cm;
248  pos_y = 0*cm;
249  pos_z = 95*cm;
250  G4VPhysicalVolume *pvol_rest =
251  new G4PVPlacement(0,
252  G4ThreeVector(pos_x, pos_y, pos_z),
253  aRest_log,
254  name,
255  worldCylinder_log,
256  false,
257  19); // i=19
258 
259  fPhysicalVolumeVector.push_back(pvol_rest);
260 
261  SetSensitive();
262  return fWorldVolume;
263 }
const XML_Char * name
Definition: expat.h:151
CLHEP::Hep3Vector G4ThreeVector
int universe_mean_density
Definition: hepunit.py:307
Definition: G4Tubs.hh:85
int G4int
Definition: G4Types.hh:78
void SetForceSolid(G4bool=true)
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
double A(double temperature)
static constexpr double cm
Definition: G4SIunits.hh:119
#define pascal
static constexpr double kelvin
Definition: G4SIunits.hh:281
static constexpr double cm3
Definition: G4SIunits.hh:121
tuple z
Definition: test.py:28
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:362
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152
static constexpr double mole
Definition: G4SIunits.hh:286
void SetVisAttributes(const G4VisAttributes *pVA)

Here is the call graph for this function:

void B01DetectorConstruction::ConstructSDandField ( )
virtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 408 of file B01DetectorConstruction.cc.

409 {
410 
411  // Sensitive Detector Manager.
413  // Sensitive Detector Name
414  G4String concreteSDname = "ConcreteSD";
415 
416  //------------------------
417  // MultiFunctionalDetector
418  //------------------------
419  //
420  // Define MultiFunctionalDetector with name.
421  G4MultiFunctionalDetector* MFDet =
422  new G4MultiFunctionalDetector(concreteSDname);
423  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
424 
425  G4String fltName,particleName;
426  G4SDParticleFilter* neutronFilter =
427  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
428 
429  MFDet->SetFilter(neutronFilter);
430 
431  for (std::vector<G4LogicalVolume *>::iterator it =
432  fLogicalVolumeVector.begin();
433  it != fLogicalVolumeVector.end(); it++){
434  // (*it)->SetSensitiveDetector(MFDet);
435  SetSensitiveDetector((*it)->GetName(), MFDet);
436  }
437 
438  G4String psName;
439  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
440  MFDet->RegisterPrimitive(scorer0);
441 
442  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
443  scorer1->Weighted(true);
444  MFDet->RegisterPrimitive(scorer1);
445 
446  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
447  MFDet->RegisterPrimitive(scorer2);
448 
449  G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter"
450  ,fCurrent_In);
451  MFDet->RegisterPrimitive(scorer3);
452 
453  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
454  MFDet->RegisterPrimitive(scorer4);
455 
456  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
457  scorer5->Weighted(true);
458  MFDet->RegisterPrimitive(scorer5);
459 
460  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
461  scorer6->Weighted(true);
462  scorer6->MultiplyKineticEnergy(true);
463  MFDet->RegisterPrimitive(scorer6);
464 
465  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
466  scorer7->Weighted(true);
467  scorer7->DivideByVelocity(true);
468  MFDet->RegisterPrimitive(scorer7);
469 
470  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
471  scorer8->Weighted(true);
472  scorer8->MultiplyKineticEnergy(true);
473  scorer8->DivideByVelocity(true);
474  MFDet->RegisterPrimitive(scorer8);
475 
476 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
void Weighted(G4bool flg=true)
void Weighted(G4bool flg=true)
void SetFilter(G4VSDFilter *value)
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void MultiplyKineticEnergy(G4bool flg=true)
void DivideByVelocity(G4bool flg=true)

Here is the call graph for this function:

G4VIStore * B01DetectorConstruction::CreateImportanceStore ( )

Definition at line 267 of file B01DetectorConstruction.cc.

268 {
269  G4cout << " B01DetectorConstruction:: Creating Importance Store " << G4endl;
270  if (!fPhysicalVolumeVector.size())
271  {
272  G4Exception("B01DetectorConstruction::CreateImportanceStore"
273  ,"exampleB01_0001",RunMustBeAborted
274  ,"no physical volumes created yet!");
275  }
276 
277  fWorldVolume = fPhysicalVolumeVector[0];
278 
279  // creating and filling the importance store
280 
281  G4IStore *istore = G4IStore::GetInstance();
282 
283  G4int n = 0;
284  G4double imp =1;
285  istore->AddImportanceGeometryCell(1, *fWorldVolume);
286  for (std::vector<G4VPhysicalVolume *>::iterator
287  it = fPhysicalVolumeVector.begin();
288  it != fPhysicalVolumeVector.end() - 1; it++)
289  {
290  if (*it != fWorldVolume)
291  {
292  imp = std::pow(2., n++);
293  G4cout << "Going to assign importance: " << imp << ", to volume: "
294  << (*it)->GetName() << G4endl;
295  istore->AddImportanceGeometryCell(imp, *(*it),n);
296  }
297  }
298 
299  // the remaining part pf the geometry (rest) gets the same
300  // importance as the last conrete cell
301  //
302  istore->AddImportanceGeometryCell(imp,
303  *(fPhysicalVolumeVector[fPhysicalVolumeVector.size()-1]),++n);
304 
305  return istore;
306 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void AddImportanceGeometryCell(G4double importance, const G4GeometryCell &gCell)
Definition: G4IStore.cc:107
const G4int n
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4IStore * GetInstance()
Definition: G4IStore.cc:236
#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:

G4VWeightWindowStore * B01DetectorConstruction::CreateWeightWindowStore ( )

Definition at line 310 of file B01DetectorConstruction.cc.

311 {
312  if (!fPhysicalVolumeVector.size())
313  {
314  G4Exception("B01DetectorConstruction::CreateWeightWindowStore"
315  ,"exampleB01_0002",RunMustBeAborted
316  ,"no physical volumes created yet!");
317  }
318 
319  fWorldVolume = fPhysicalVolumeVector[0];
320 
321  // creating and filling the weight window store
322 
324 
325  // create one energy region covering the energies of the problem
326  //
327  std::set<G4double, std::less<G4double> > enBounds;
328  enBounds.insert(1 * GeV);
329  wwstore->SetGeneralUpperEnergyBounds(enBounds);
330 
331  G4int n = 0;
332  G4double lowerWeight =1;
333  std::vector<G4double> lowerWeights;
334 
335  lowerWeights.push_back(1);
336  G4GeometryCell gWorldCell(*fWorldVolume,0);
337  wwstore->AddLowerWeights(gWorldCell, lowerWeights);
338 
339  for (std::vector<G4VPhysicalVolume *>::iterator
340  it = fPhysicalVolumeVector.begin();
341  it != fPhysicalVolumeVector.end() - 1; it++)
342  {
343  if (*it != fWorldVolume)
344  {
345  lowerWeight = 1./std::pow(2., n++);
346  G4cout << "Going to assign lower weight: " << lowerWeight
347  << ", to volume: "
348  << (*it)->GetName() << G4endl;
349  G4GeometryCell gCell(*(*it),n);
350  lowerWeights.clear();
351  lowerWeights.push_back(lowerWeight);
352  wwstore->AddLowerWeights(gCell, lowerWeights);
353  }
354  }
355 
356  // the remaining part pf the geometry (rest) gets the same
357  // lower weight bound as the last conrete cell
358  //
360  gRestCell(*(fPhysicalVolumeVector[fPhysicalVolumeVector.size()-1]), ++n);
361  wwstore->AddLowerWeights(gRestCell, lowerWeights);
362 
363  return wwstore;
364 }
void AddLowerWeights(const G4GeometryCell &gCell, const std::vector< G4double > &lowerWeights)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4int n
static G4WeightWindowStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetGeneralUpperEnergyBounds(const std::set< G4double, std::less< G4double > > &enBounds)
static constexpr double GeV
Definition: G4SIunits.hh:217
#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:

G4String B01DetectorConstruction::GetCellName ( G4int  i)

Definition at line 368 of file B01DetectorConstruction.cc.

369 {
370  std::ostringstream os;
371  os << "cell_";
372  if (i<10)
373  {
374  os << "0";
375  }
376  os << i ;
377  G4String name = os.str();
378  return name;
379 }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

G4VPhysicalVolume * B01DetectorConstruction::GetWorldVolume ( )

Definition at line 381 of file B01DetectorConstruction.cc.

381  {
382  return fWorldVolume;
383 }

Here is the caller graph for this function:

void B01DetectorConstruction::SetSensitive ( )

Definition at line 387 of file B01DetectorConstruction.cc.

387  {
388 
389  // -------------------------------------------------
390  // The collection names of defined Primitives are
391  // 0 ConcreteSD/Collisions
392  // 1 ConcreteSD/CollWeight
393  // 2 ConcreteSD/Population
394  // 3 ConcreteSD/TrackEnter
395  // 4 ConcreteSD/SL
396  // 5 ConcreteSD/SLW
397  // 6 ConcreteSD/SLWE
398  // 7 ConcreteSD/SLW_V
399  // 8 ConcreteSD/SLWE_V
400  // -------------------------------------------------
401 
402  // moved to ConstructSDandField() for MT compliance
403 
404 }

Here is the caller graph for this function:


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