Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B02ImportanceDetectorConstruction.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 "globals.hh"
34 #include <sstream>
35 
37 
38 #include "G4Material.hh"
39 #include "G4Tubs.hh"
40 #include "G4LogicalVolume.hh"
41 #include "G4ThreeVector.hh"
42 #include "G4PVPlacement.hh"
43 #include "G4PhysicalConstants.hh"
44 #include "G4SystemOfUnits.hh"
45 
46 // For Primitive Scorers
47 #include "G4SDManager.hh"
49 #include "G4SDParticleFilter.hh"
50 #include "G4PSNofCollision.hh"
51 #include "G4PSPopulation.hh"
52 #include "G4PSTrackCounter.hh"
53 #include "G4PSTrackLength.hh"
54 
55 
56 
58 :G4VUserParallelWorld(worldName),fLogicalVolumeVector()
59 {
60  // Construct();
61 }
62 
64 {
65  fLogicalVolumeVector.clear();
66 }
67 
68 void B02ImportanceDetectorConstruction::Construct()
69 {
70  G4cout << " constructing parallel world " << G4endl;
71 
72  //GetWorld methods create a clone of the mass world to the parallel world (!)
73  // via the transportation manager
74  ghostWorld = GetWorld();
75  G4LogicalVolume* worldLogical = ghostWorld->GetLogicalVolume();
76  fLogicalVolumeVector.push_back(worldLogical);
77 
78  G4String name("none");
79  G4double density(universe_mean_density), temperature(0), pressure(0);
80 
81  name = "Galactic";
82  density = universe_mean_density; //from PhysicalConstants.h
83  pressure = 3.e-18*pascal;
84  temperature = 2.73*kelvin;
85  G4cout << density << " " << kStateGas << G4endl;
86  G4Material *Galactic =
87  new G4Material(name, 1., 1.01*g/mole, density,
88  kStateGas,temperature,pressure);
89 
90 
91  // fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
92  fPVolumeStore.AddPVolume(G4GeometryCell(*ghostWorld, 0));
93 
94 
95 
96 
97  // creating 18 slobs of 10 cm thicknes
98 
99  G4double innerRadiusShield = 0*cm;
100  G4double outerRadiusShield = 100*cm;
101  G4double hightShield = 5*cm;
102  G4double startAngleShield = 0*deg;
103  G4double spanningAngleShield = 360*deg;
104 
105  G4Tubs *aShield = new G4Tubs("aShield",
106  innerRadiusShield,
107  outerRadiusShield,
108  hightShield,
109  startAngleShield,
110  spanningAngleShield);
111 
112  // logical parallel cells
113 
114  G4LogicalVolume *aShield_log =
115  new G4LogicalVolume(aShield, Galactic, "aShield_log");
116  fLogicalVolumeVector.push_back(aShield_log);
117 
118  // physical parallel cells
119 
120  G4int i = 1;
121  G4double startz = -85*cm;
122  // for (i=1; i<=18; ++i) {
123  for (i=1; i<=18; i++) {
124 
125  name = GetCellName(i);
126 
127  G4double pos_x = 0*cm;
128  G4double pos_y = 0*cm;
129  G4double pos_z = startz + (i-1) * (2*hightShield);
130  G4VPhysicalVolume *pvol =
131  new G4PVPlacement(0,
132  G4ThreeVector(pos_x, pos_y, pos_z),
133  aShield_log,
134  name,
135  worldLogical,
136  false,
137  i);
138  // 0);
139  G4GeometryCell cell(*pvol, i);
140  // G4GeometryCell cell(*pvol, 0);
141  fPVolumeStore.AddPVolume(cell);
142  }
143 
144  // filling the rest of the world volumr behind the concrete with
145  // another slob which should get the same importance value as the
146  // last slob
147  innerRadiusShield = 0*cm;
148  // outerRadiusShield = 110*cm; exceeds world volume!!!!
149  outerRadiusShield = 101*cm;
150  // hightShield = 10*cm;
151  hightShield = 5*cm;
152  startAngleShield = 0*deg;
153  spanningAngleShield = 360*deg;
154 
155  G4Tubs *aRest = new G4Tubs("Rest",
156  innerRadiusShield,
157  outerRadiusShield,
158  hightShield,
159  startAngleShield,
160  spanningAngleShield);
161 
162  G4LogicalVolume *aRest_log =
163  new G4LogicalVolume(aRest, Galactic, "aRest_log");
164 
165  fLogicalVolumeVector.push_back(aRest_log);
166 
167  name = GetCellName(19);
168 
169  G4double pos_x = 0*cm;
170  G4double pos_y = 0*cm;
171  // G4double pos_z = 100*cm;
172  G4double pos_z = 95*cm;
173  G4VPhysicalVolume *pvol =
174  new G4PVPlacement(0,
175  G4ThreeVector(pos_x, pos_y, pos_z),
176  aRest_log,
177  name,
178  worldLogical,
179  false,
180  19);
181  // 0);
182  G4GeometryCell cell(*pvol, 19);
183  // G4GeometryCell cell(*pvol, 0);
184  fPVolumeStore.AddPVolume(cell);
185 
186  SetSensitive();
187 
188 }
189 
192  return *fPVolumeStore.GetPVolume(name);
193 }
194 
195 
197  G4String names(fPVolumeStore.GetPNames());
198  return names;
199 }
200 
201 
203  std::ostringstream os;
204  os << "cell_";
205  if (i<10) {
206  os << "0";
207  }
208  os << i;
209  G4String name = os.str();
210  return name;
211 }
212 
215  const G4VPhysicalVolume *p=0;
216  p = fPVolumeStore.GetPVolume(name);
217  if (p) {
218  return G4GeometryCell(*p,0);
219  }
220  else {
221  G4cout << "B02ImportanceDetectorConstruction::GetGeometryCell: couldn't get G4GeometryCell" << G4endl;
222  return G4GeometryCell(*ghostWorld,-2);
223  }
224 }
225 
226 
228  return *ghostWorld;
229 }
230 
232  return ghostWorld;
233 }
234 
235 
237 
238  // -------------------------------------------------
239  // The collection names of defined Primitives are
240  // 0 ConcreteSD/Collisions
241  // 1 ConcreteSD/CollWeight
242  // 2 ConcreteSD/Population
243  // 3 ConcreteSD/TrackEnter
244  // 4 ConcreteSD/SL
245  // 5 ConcreteSD/SLW
246  // 6 ConcreteSD/SLWE
247  // 7 ConcreteSD/SLW_V
248  // 8 ConcreteSD/SLWE_V
249  // -------------------------------------------------
250 
251 
252  //================================================
253  // Sensitive detectors : MultiFunctionalDetector
254  //================================================
255  //
256  // Sensitive Detector Manager.
257 
259  //
260  // Sensitive Detector Name
261  G4String concreteSDname = "ConcreteSD";
262 
263  //------------------------
264  // MultiFunctionalDetector
265  //------------------------
266  //
267  // Define MultiFunctionalDetector with name.
268  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
269  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
270 
271 
272  G4String fltName,particleName;
273  G4SDParticleFilter* neutronFilter =
274  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
275 
276  MFDet->SetFilter(neutronFilter);
277 
278 
279  for (std::vector<G4LogicalVolume *>::iterator it = fLogicalVolumeVector.begin();
280  it != fLogicalVolumeVector.end(); it++){
281  (*it)->SetSensitiveDetector(MFDet);
282  }
283 
284  G4String psName;
285  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
286  MFDet->RegisterPrimitive(scorer0);
287 
288 
289  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
290  scorer1->Weighted(true);
291  MFDet->RegisterPrimitive(scorer1);
292 
293 
294  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
295  MFDet->RegisterPrimitive(scorer2);
296 
297  G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
298  MFDet->RegisterPrimitive(scorer3);
299 
300  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
301  MFDet->RegisterPrimitive(scorer4);
302 
303  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
304  scorer5->Weighted(true);
305  MFDet->RegisterPrimitive(scorer5);
306 
307  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
308  scorer6->Weighted(true);
309  scorer6->MultiplyKineticEnergy(true);
310  MFDet->RegisterPrimitive(scorer6);
311 
312  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
313  scorer7->Weighted(true);
314  scorer7->DivideByVelocity(true);
315  MFDet->RegisterPrimitive(scorer7);
316 
317  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
318  scorer8->Weighted(true);
319  scorer8->MultiplyKineticEnergy(true);
320  scorer8->DivideByVelocity(true);
321  MFDet->RegisterPrimitive(scorer8);
322 
323 }