Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B02PSScoringDetectorConstruction.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 B02PSScoringDetectorConstruction::Construct()
69 {
70 
71  G4cout << " constructing parallel world " << G4endl;
72 
73  //GetWorld methods create a clone of the mass world to the parallel world (!)
74  // via the transportation manager
75  ghostWorld = GetWorld();
76  G4LogicalVolume* worldLogical = ghostWorld->GetLogicalVolume();
77  fLogicalVolumeVector.push_back(worldLogical);
78 
79  G4String name("none");
80  G4double density(universe_mean_density), temperature(0), pressure(0);
81 
82  name = "Galactic";
83  density = universe_mean_density; //from PhysicalConstants.h
84  pressure = 3.e-18*pascal;
85  temperature = 2.73*kelvin;
86  G4cout << density << " " << kStateGas << G4endl;
87  G4Material *Galactic =
88  new G4Material(name, 1., 1.01*g/mole, density,
89  kStateGas,temperature,pressure);
90 
91 
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 
124  name = GetCellName(i);
125 
126  G4double pos_x = 0*cm;
127  G4double pos_y = 0*cm;
128  G4double pos_z = startz + (i-1) * (2*hightShield);
129  new G4PVPlacement(0,
130  G4ThreeVector(pos_x, pos_y, pos_z),
131  aShield_log,
132  name,
133  worldLogical,
134  false,
135  i);
136  //G4GeometryCell cell(*pvol, i);
137  //fPVolumeStore.AddPVolume(cell);
138  }
139 
140  // filling the rest of the world volumr behind the concrete with
141  // another slob which should get the same importance value as the
142  // last slob
143  innerRadiusShield = 0*cm;
144  outerRadiusShield = 100*cm;
145  hightShield = 5*cm;
146  startAngleShield = 0*deg;
147  spanningAngleShield = 360*deg;
148 
149  G4Tubs *aRest = new G4Tubs("Rest",
150  innerRadiusShield,
151  outerRadiusShield,
152  hightShield,
153  startAngleShield,
154  spanningAngleShield);
155 
156  G4LogicalVolume *aRest_log =
157  new G4LogicalVolume(aRest, Galactic, "aRest_log");
158  name = GetCellName(19);
159  fLogicalVolumeVector.push_back(aRest_log);
160 
161  G4double pos_x = 0*cm;
162  G4double pos_y = 0*cm;
163  G4double pos_z = 95*cm;
164  new G4PVPlacement(0,
165  G4ThreeVector(pos_x, pos_y, pos_z),
166  aRest_log,
167  name,
168  worldLogical,
169  false,
170  19); // i = 19
171  //G4GeometryCell cell(*pvol, i);
172  //fPVolumeStore.AddPVolume(cell);
173 
174  SetSensitive();
175 
176 }
177 
179  std::ostringstream os;
180  os << "cell_";
181  if (i<10) {
182  os << "0";
183  }
184  os << i;
185  G4String name = os.str();
186  return name;
187 }
188 
190  return ghostWorld;
191 }
192 
193 
195  return *ghostWorld;
196 }
197 
199 
200  // -------------------------------------------------
201  // The collection names of defined Primitives are
202  // 0 ConcreteSD/Collisions
203  // 1 ConcreteSD/CollWeight
204  // 2 ConcreteSD/Population
205  // 3 ConcreteSD/TrackEnter
206  // 4 ConcreteSD/SL
207  // 5 ConcreteSD/SLW
208  // 6 ConcreteSD/SLWE
209  // 7 ConcreteSD/SLW_V
210  // 8 ConcreteSD/SLWE_V
211  // -------------------------------------------------
212 
213 
214  //================================================
215  // Sensitive detectors : MultiFunctionalDetector
216  //================================================
217  //
218  // Sensitive Detector Manager.
220  //
221  // Sensitive Detector Name
222  G4String concreteSDname = "ConcreteSD";
223 
224  //------------------------
225  // MultiFunctionalDetector
226  //------------------------
227  //
228  // Define MultiFunctionalDetector with name.
229  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
230  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
231 
232 
233  G4String fltName,particleName;
234  G4SDParticleFilter* neutronFilter =
235  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
236 
237  MFDet->SetFilter(neutronFilter);
238 
239 
240  for (std::vector<G4LogicalVolume *>::iterator it = fLogicalVolumeVector.begin();
241  it != fLogicalVolumeVector.end(); it++){
242  (*it)->SetSensitiveDetector(MFDet);
243  }
244 
245  G4String psName;
246  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
247  MFDet->RegisterPrimitive(scorer0);
248 
249 
250  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
251  scorer1->Weighted(true);
252  MFDet->RegisterPrimitive(scorer1);
253 
254 
255  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
256  MFDet->RegisterPrimitive(scorer2);
257 
258  G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
259  MFDet->RegisterPrimitive(scorer3);
260 
261  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
262  MFDet->RegisterPrimitive(scorer4);
263 
264  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
265  scorer5->Weighted(true);
266  MFDet->RegisterPrimitive(scorer5);
267 
268  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
269  scorer6->Weighted(true);
270  scorer6->MultiplyKineticEnergy(true);
271  MFDet->RegisterPrimitive(scorer6);
272 
273  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
274  scorer7->Weighted(true);
275  scorer7->DivideByVelocity(true);
276  MFDet->RegisterPrimitive(scorer7);
277 
278  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
279  scorer8->Weighted(true);
280  scorer8->MultiplyKineticEnergy(true);
281  scorer8->DivideByVelocity(true);
282  MFDet->RegisterPrimitive(scorer8);
283 
284 }