Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B01DetectorConstruction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
28 //
29 //
30 // $Id$
31 //
32 
33 #include "G4Types.hh"
34 #include <sstream>
35 #include <set>
36 #include "globals.hh"
37 
39 
40 #include "G4Material.hh"
41 #include "G4Box.hh"
42 #include "G4Tubs.hh"
43 #include "G4LogicalVolume.hh"
44 #include "G4ThreeVector.hh"
45 #include "G4PVPlacement.hh"
46 #include "G4VisAttributes.hh"
47 #include "G4Colour.hh"
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 
51 // For Primitive Scorers
52 #include "G4SDManager.hh"
54 #include "G4SDParticleFilter.hh"
55 #include "G4PSNofCollision.hh"
56 #include "G4PSPopulation.hh"
57 #include "G4PSTrackCounter.hh"
58 #include "G4PSTrackLength.hh"
59 
60 // for importance biasing
61 #include "G4IStore.hh"
62 
63 // for weight window technique
64 #include "G4WeightWindowStore.hh"
65 
67  fPhysicalVolumeVector(),fLogicalVolumeVector()
68 {;}
69 
71 {;}
72 
74 {
75  G4double pos_x;
76  G4double pos_y;
77  G4double pos_z;
78 
79  G4double density, pressure, temperature;
80  G4double A;
81  G4int Z;
82 
84  G4double z;
85  G4double fractionmass;
86 
87  A = 1.01*g/mole;
88  G4Element* elH = new G4Element(name="Hydrogen",symbol="H" , Z= 1, A);
89 
90  A = 12.01*g/mole;
91  G4Element* elC = new G4Element(name="Carbon" ,symbol="C" , Z = 6, A);
92 
93  A = 16.00*g/mole;
94  G4Element* elO = new G4Element(name="Oxygen" ,symbol="O" , Z= 8, A);
95 
96  A = 22.99*g/mole;
97  G4Element* elNa = new G4Element(name="Natrium" ,symbol="Na" , Z=11 , A);
98 
99  A = 200.59*g/mole;
100  G4Element* elHg = new G4Element(name="Hg" ,symbol="Hg" , Z=80, A);
101 
102  A = 26.98*g/mole;
103  G4Element* elAl = new G4Element(name="Aluminium" ,symbol="Al" , Z=13, A);
104 
105  A = 28.09*g/mole;
106  G4Element* elSi = new G4Element(name="Silicon", symbol="Si", Z=14, A);
107 
108  A = 39.1*g/mole;
109  G4Element* elK = new G4Element(name="K" ,symbol="K" , Z=19 , A);
110 
111  A = 69.72*g/mole;
112  G4Element* elCa = new G4Element(name="Calzium" ,symbol="Ca" , Z=31 , A);
113 
114  A = 55.85*g/mole;
115  G4Element* elFe = new G4Element(name="Iron" ,symbol="Fe", Z=26, A);
116 
117  density = universe_mean_density; //from PhysicalConstants.h
118  pressure = 3.e-18*pascal;
119  temperature = 2.73*kelvin;
120  G4Material *Galactic =
121  new G4Material(name="Galactic", z=1., A=1.01*g/mole, density,
122  kStateGas,temperature,pressure);
123 
124  density = 2.03*g/cm3;
125  G4Material* Concrete = new G4Material("Concrete", density, 10);
126  Concrete->AddElement(elH , fractionmass= 0.01);
127  Concrete->AddElement(elO , fractionmass= 0.529);
128  Concrete->AddElement(elNa , fractionmass= 0.016);
129  Concrete->AddElement(elHg , fractionmass= 0.002);
130  Concrete->AddElement(elAl , fractionmass= 0.034);
131  Concrete->AddElement(elSi , fractionmass= 0.337);
132  Concrete->AddElement(elK , fractionmass= 0.013);
133  Concrete->AddElement(elCa , fractionmass= 0.044);
134  Concrete->AddElement(elFe , fractionmass= 0.014);
135  Concrete->AddElement(elC , fractionmass= 0.001);
136 
138  // world cylinder volume
140 
141  // world solid
142 
143  G4double innerRadiusCylinder = 0*cm;
144  G4double outerRadiusCylinder = 100*cm;
145  G4double hightCylinder = 100*cm;
146  G4double startAngleCylinder = 0*deg;
147  G4double spanningAngleCylinder = 360*deg;
148 
149  G4Tubs *worldCylinder = new G4Tubs("worldCylinder",
150  innerRadiusCylinder,
151  outerRadiusCylinder,
152  hightCylinder,
153  startAngleCylinder,
154  spanningAngleCylinder);
155 
156  // logical world
157 
158  G4LogicalVolume *worldCylinder_log =
159  new G4LogicalVolume(worldCylinder, Galactic, "worldCylinder_log");
160  fLogicalVolumeVector.push_back(worldCylinder_log);
161 
162  name = "shieldWorld";
163  pWorldVolume = new
164  G4PVPlacement(0, G4ThreeVector(0,0,0), worldCylinder_log,
165  name, 0, false, 0);
166 
167 
168  fPhysicalVolumeVector.push_back(pWorldVolume);
169 
170  // creating 18 slabs of 10 cm thick concrete
171 
172  G4double innerRadiusShield = 0*cm;
173  G4double outerRadiusShield = 100*cm;
174  G4double hightShield = 5*cm;
175  G4double startAngleShield = 0*deg;
176  G4double spanningAngleShield = 360*deg;
177 
178  G4Tubs *aShield = new G4Tubs("aShield",
179  innerRadiusShield,
180  outerRadiusShield,
181  hightShield,
182  startAngleShield,
183  spanningAngleShield);
184 
185  // logical shield
186 
187  G4LogicalVolume *aShield_log =
188  new G4LogicalVolume(aShield, Concrete, "aShield_log");
189  fLogicalVolumeVector.push_back(aShield_log);
190 
191  G4VisAttributes* pShieldVis = new G4VisAttributes(G4Colour(0.0,0.0,1.0));
192  pShieldVis->SetForceSolid(true);
193  aShield_log->SetVisAttributes(pShieldVis);
194 
195  // physical shields
196 
197  G4int i;
198  G4double startz = -85*cm;
199  for (i=1; i<=18; i++)
200  {
201  name = GetCellName(i);
202  pos_x = 0*cm;
203  pos_y = 0*cm;
204  pos_z = startz + (i-1) * (2*hightShield);
205  G4VPhysicalVolume *pvol =
206  new G4PVPlacement(0,
207  G4ThreeVector(pos_x, pos_y, pos_z),
208  aShield_log,
209  name,
210  worldCylinder_log,
211  false,
212  i);
213  fPhysicalVolumeVector.push_back(pvol);
214  }
215 
216  // filling the rest of the world volume behind the concrete with
217  // another slab which should get the same importance value
218  // or lower weight bound as the last slab
219  //
220  innerRadiusShield = 0*cm;
221  outerRadiusShield = 100*cm;
222  hightShield = 5*cm;
223  startAngleShield = 0*deg;
224  spanningAngleShield = 360*deg;
225 
226  G4Tubs *aRest = new G4Tubs("Rest",
227  innerRadiusShield,
228  outerRadiusShield,
229  hightShield,
230  startAngleShield,
231  spanningAngleShield);
232 
233  G4LogicalVolume *aRest_log =
234  new G4LogicalVolume(aRest, Galactic, "aRest_log");
235  fLogicalVolumeVector.push_back(aRest_log);
236  name = "rest";
237 
238  pos_x = 0*cm;
239  pos_y = 0*cm;
240  pos_z = 95*cm;
241  G4VPhysicalVolume *pvol_rest =
242  new G4PVPlacement(0,
243  G4ThreeVector(pos_x, pos_y, pos_z),
244  aRest_log,
245  name,
246  worldCylinder_log,
247  false,
248  19); // i=19
249 
250  fPhysicalVolumeVector.push_back(pvol_rest);
251 
252  SetSensitive();
253  return pWorldVolume;
254 }
255 
256 
258 {
259  if (!fPhysicalVolumeVector.size())
260  {
261  G4Exception("B01DetectorConstruction::CreateImportanceStore","exampleB01_0001",RunMustBeAborted,"no physical volumes created yet!");
262  }
263 
264  pWorldVolume = fPhysicalVolumeVector[0];
265 
266  // creating and filling the importance store
267 
268  G4IStore *istore = new G4IStore(*pWorldVolume);
269 
270  G4int n = 0;
271  G4double imp =1;
272  istore->AddImportanceGeometryCell(1, *pWorldVolume);
273  for (std::vector<G4VPhysicalVolume *>::iterator
274  it = fPhysicalVolumeVector.begin();
275  it != fPhysicalVolumeVector.end() - 1; it++)
276  {
277  if (*it != pWorldVolume)
278  {
279  imp = std::pow(2., n++);
280  G4cout << "Going to assign importance: " << imp << ", to volume: "
281  << (*it)->GetName() << G4endl;
282  istore->AddImportanceGeometryCell(imp, *(*it),n);
283  }
284  }
285 
286  // the remaining part pf the geometry (rest) gets the same
287  // importance as the last conrete cell
288  //
289  istore->AddImportanceGeometryCell(imp,
290  *(fPhysicalVolumeVector[fPhysicalVolumeVector.size()-1]),++n);
291 
292  return istore;
293 }
294 
296 {
297  if (!fPhysicalVolumeVector.size())
298  {
299  G4Exception("B01DetectorConstruction::CreateWeightWindowStore","exampleB01_0002",RunMustBeAborted,"no physical volumes created yet!");
300  }
301 
302  pWorldVolume = fPhysicalVolumeVector[0];
303 
304  // creating and filling the weight window store
305 
306  G4WeightWindowStore *wwstore = new G4WeightWindowStore(*pWorldVolume);
307 
308  // create one energy region covering the energies of the problem
309  //
310  std::set<G4double, std::less<G4double> > enBounds;
311  enBounds.insert(1 * GeV);
312  wwstore->SetGeneralUpperEnergyBounds(enBounds);
313 
314  G4int n = 0;
315  G4double lowerWeight =1;
316  std::vector<G4double> lowerWeights;
317 
318  lowerWeights.push_back(1);
319  G4GeometryCell gWorldCell(*pWorldVolume,0);
320  wwstore->AddLowerWeights(gWorldCell, lowerWeights);
321 
322  for (std::vector<G4VPhysicalVolume *>::iterator
323  it = fPhysicalVolumeVector.begin();
324  it != fPhysicalVolumeVector.end() - 1; it++)
325  {
326  if (*it != pWorldVolume)
327  {
328  lowerWeight = 1./std::pow(2., n++);
329  G4cout << "Going to assign lower weight: " << lowerWeight
330  << ", to volume: "
331  << (*it)->GetName() << G4endl;
332  G4GeometryCell gCell(*(*it),n);
333  lowerWeights.clear();
334  lowerWeights.push_back(lowerWeight);
335  wwstore->AddLowerWeights(gCell, lowerWeights);
336  }
337  }
338 
339  // the remaining part pf the geometry (rest) gets the same
340  // lower weight bound as the last conrete cell
341  //
343  gRestCell(*(fPhysicalVolumeVector[fPhysicalVolumeVector.size()-1]), ++n);
344  wwstore->AddLowerWeights(gRestCell, lowerWeights);
345 
346  return wwstore;
347 }
348 
350 {
351  std::ostringstream os;
352  os << "cell_";
353  if (i<10)
354  {
355  os << "0";
356  }
357  os << i ;
358  G4String name = os.str();
359  return name;
360 }
361 
363  return pWorldVolume;
364 }
365 
367 
368  // -------------------------------------------------
369  // The collection names of defined Primitives are
370  // 0 ConcreteSD/Collisions
371  // 1 ConcreteSD/CollWeight
372  // 2 ConcreteSD/Population
373  // 3 ConcreteSD/TrackEnter
374  // 4 ConcreteSD/SL
375  // 5 ConcreteSD/SLW
376  // 6 ConcreteSD/SLWE
377  // 7 ConcreteSD/SLW_V
378  // 8 ConcreteSD/SLWE_V
379  // -------------------------------------------------
380 
381 
382  //================================================
383  // Sensitive detectors : MultiFunctionalDetector
384  //================================================
385  //
386  // Sensitive Detector Manager.
388  //
389  // Sensitive Detector Name
390  G4String concreteSDname = "ConcreteSD";
391 
392  //------------------------
393  // MultiFunctionalDetector
394  //------------------------
395  //
396  // Define MultiFunctionalDetector with name.
397  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
398  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
399 
400 
401  G4String fltName,particleName;
402  G4SDParticleFilter* neutronFilter =
403  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
404 
405  MFDet->SetFilter(neutronFilter);
406 
407 
408  for (std::vector<G4LogicalVolume *>::iterator it = fLogicalVolumeVector.begin();
409  it != fLogicalVolumeVector.end(); it++){
410  (*it)->SetSensitiveDetector(MFDet);
411  }
412 
413  G4String psName;
414  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
415  MFDet->RegisterPrimitive(scorer0);
416 
417 
418  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
419  scorer1->Weighted(true);
420  MFDet->RegisterPrimitive(scorer1);
421 
422 
423  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
424  MFDet->RegisterPrimitive(scorer2);
425 
426  G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
427  MFDet->RegisterPrimitive(scorer3);
428 
429  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
430  MFDet->RegisterPrimitive(scorer4);
431 
432  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
433  scorer5->Weighted(true);
434  MFDet->RegisterPrimitive(scorer5);
435 
436  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
437  scorer6->Weighted(true);
438  scorer6->MultiplyKineticEnergy(true);
439  MFDet->RegisterPrimitive(scorer6);
440 
441  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
442  scorer7->Weighted(true);
443  scorer7->DivideByVelocity(true);
444  MFDet->RegisterPrimitive(scorer7);
445 
446  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
447  scorer8->Weighted(true);
448  scorer8->MultiplyKineticEnergy(true);
449  scorer8->DivideByVelocity(true);
450  MFDet->RegisterPrimitive(scorer8);
451 
452 }