Geant4_10
B03ImportanceDetectorConstruction.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: B03ImportanceDetectorConstruction.cc 70093 2013-05-23 09:03:57Z gcosmo $
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
60 :G4VUserParallelWorld(worldName),fLogicalVolumeVector()
61 {
62  // Construct();
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
68 {
69  fLogicalVolumeVector.clear();
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
74 void B03ImportanceDetectorConstruction::Construct()
75 {
76  G4cout << " constructing parallel world " << G4endl;
77 
78  G4Material* dummyMat = 0;
79 
80  //GetWorld methods create a clone of the mass world to the parallel world (!)
81  // via the transportation manager
82  fGhostWorld = GetWorld();
83  G4cout << " B03ImportanceDetectorConstruction:: ghostWorldName = "
84  << fGhostWorld->GetName() << G4endl;
85  G4LogicalVolume* worldLogical = fGhostWorld->GetLogicalVolume();
86  fLogicalVolumeVector.push_back(worldLogical);
87 
88  G4String name("none");
89  // fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
90  fPVolumeStore.AddPVolume(G4GeometryCell(*fGhostWorld, 0));
91 
92 
93 
94 
95  // creating 18 slobs of 10 cm thicknes
96 
97  G4double innerRadiusShield = 0*cm;
98  G4double outerRadiusShield = 100*cm;
99  G4double heightShield = 5*cm;
100  G4double startAngleShield = 0*deg;
101  G4double spanningAngleShield = 360*deg;
102 
103  G4Tubs *aShield = new G4Tubs("aShield",
104  innerRadiusShield,
105  outerRadiusShield,
106  heightShield,
107  startAngleShield,
108  spanningAngleShield);
109 
110  // logical parallel cells
111 
112  G4LogicalVolume *aShield_log_imp =
113  new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
114  fLogicalVolumeVector.push_back(aShield_log_imp);
115 
116  // physical parallel cells
117 
118  G4int i = 1;
119  G4double startz = -85*cm;
120  // for (i=1; i<=18; ++i) {
121  for (i=1; i<=18; i++) {
122 
123  name = GetCellName(i);
124 
125  G4double pos_x = 0*cm;
126  G4double pos_y = 0*cm;
127  G4double pos_z = startz + (i-1) * (2*heightShield);
128  G4VPhysicalVolume *pvol =
129  new G4PVPlacement(0,
130  G4ThreeVector(pos_x, pos_y, pos_z),
131  aShield_log_imp,
132  name,
133  worldLogical,
134  false,
135  i);
136  // 0);
137  G4GeometryCell cell(*pvol, i);
138  // G4GeometryCell cell(*pvol, 0);
139  fPVolumeStore.AddPVolume(cell);
140  }
141 
142  // filling the rest of the world volumr behind the concrete with
143  // another slob which should get the same importance value as the
144  // last slob
145  innerRadiusShield = 0*cm;
146  // outerRadiusShield = 110*cm; exceeds world volume!!!!
147  outerRadiusShield = 100*cm;
148  // heightShield = 10*cm;
149  heightShield = 5*cm;
150  startAngleShield = 0*deg;
151  spanningAngleShield = 360*deg;
152 
153  G4Tubs *aRest = new G4Tubs("Rest",
154  innerRadiusShield,
155  outerRadiusShield,
156  heightShield,
157  startAngleShield,
158  spanningAngleShield);
159 
160  G4LogicalVolume *aRest_log =
161  new G4LogicalVolume(aRest, dummyMat, "aRest_log");
162 
163  fLogicalVolumeVector.push_back(aRest_log);
164 
165  name = GetCellName(19);
166 
167  G4double pos_x = 0*cm;
168  G4double pos_y = 0*cm;
169  // G4double pos_z = 100*cm;
170  G4double pos_z = 95*cm;
171  G4VPhysicalVolume *pvol =
172  new G4PVPlacement(0,
173  G4ThreeVector(pos_x, pos_y, pos_z),
174  aRest_log,
175  name,
176  worldLogical,
177  false,
178  19);
179  // 0);
180  G4GeometryCell cell(*pvol, 19);
181  // G4GeometryCell cell(*pvol, 0);
182  fPVolumeStore.AddPVolume(cell);
183 
184  SetSensitive();
185 
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
192  return *fPVolumeStore.GetPVolume(name);
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
198  G4String names(fPVolumeStore.GetPNames());
199  return names;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203 
205  std::ostringstream os;
206  os << "cell_";
207  if (i<10) {
208  os << "0";
209  }
210  os << i;
211  G4String name = os.str();
212  return name;
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 
219  const G4VPhysicalVolume *p=0;
220  p = fPVolumeStore.GetPVolume(name);
221  if (p) {
222  return G4GeometryCell(*p,0);
223  }
224  else {
225  G4cout << "B03ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
226  << " couldn't get G4GeometryCell" << G4endl;
227  return G4GeometryCell(*fGhostWorld,-2);
228  }
229 }
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
232 
235  return *fGhostWorld;
236 }
237 
238 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
239 
241  return fGhostWorld;
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245 
247 
248  // -------------------------------------------------
249  // The collection names of defined Primitives are
250  // 0 ConcreteSD/Collisions
251  // 1 ConcreteSD/CollWeight
252  // 2 ConcreteSD/Population
253  // 3 ConcreteSD/TrackEnter
254  // 4 ConcreteSD/SL
255  // 5 ConcreteSD/SLW
256  // 6 ConcreteSD/SLWE
257  // 7 ConcreteSD/SLW_V
258  // 8 ConcreteSD/SLWE_V
259  // -------------------------------------------------
260 
261  //now moved to ConstructSD()
262 
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267 {
268 
270  //
271  // Sensitive Detector Name
272  G4String concreteSDname = "ConcreteSD";
273 
274  //------------------------
275  // MultiFunctionalDetector
276  //------------------------
277  //
278  // Define MultiFunctionalDetector with name.
279  G4MultiFunctionalDetector* MFDet =
280  new G4MultiFunctionalDetector(concreteSDname);
281  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
282 
283 
284  G4String fltName,particleName;
285  G4SDParticleFilter* neutronFilter =
286  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
287 
288  MFDet->SetFilter(neutronFilter);
289 
290 
291  for (std::vector<G4LogicalVolume *>::iterator it =
292  fLogicalVolumeVector.begin();
293  it != fLogicalVolumeVector.end(); it++){
294  // (*it)->SetSensitiveDetector(MFDet);
295  SetSensitiveDetector((*it)->GetName(), MFDet);
296  }
297 
298  G4String psName;
299  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
300  MFDet->RegisterPrimitive(scorer0);
301 
302 
303  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
304  scorer1->Weighted(true);
305  MFDet->RegisterPrimitive(scorer1);
306 
307 
308  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
309  MFDet->RegisterPrimitive(scorer2);
310 
311  G4PSTrackCounter* scorer3 =
312  new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
313  MFDet->RegisterPrimitive(scorer3);
314 
315  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
316  MFDet->RegisterPrimitive(scorer4);
317 
318  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
319  scorer5->Weighted(true);
320  MFDet->RegisterPrimitive(scorer5);
321 
322  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
323  scorer6->Weighted(true);
324  scorer6->MultiplyKineticEnergy(true);
325  MFDet->RegisterPrimitive(scorer6);
326 
327  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
328  scorer7->Weighted(true);
329  scorer7->DivideByVelocity(true);
330  MFDet->RegisterPrimitive(scorer7);
331 
332  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
333  scorer8->Weighted(true);
334  scorer8->MultiplyKineticEnergy(true);
335  scorer8->DivideByVelocity(true);
336  MFDet->RegisterPrimitive(scorer8);
337 }
338 
339 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
CLHEP::Hep3Vector G4ThreeVector
G4VPhysicalVolume * GetWorld()
const char * p
Definition: xmltok.h:285
Definition: G4Tubs.hh:84
const XML_Char * name
Definition: expat.h:151
int G4int
Definition: G4Types.hh:78
void Weighted(G4bool flg=true)
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
void Weighted(G4bool flg=true)
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
void SetFilter(G4VSDFilter *value)
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:67
G4LogicalVolume * GetLogicalVolume() const
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
const G4VPhysicalVolume & GetPhysicalVolumeByName(const G4String &name) const
void MultiplyKineticEnergy(G4bool flg=true)
#define G4endl
Definition: G4ios.hh:61
G4String GetPNames() const
double G4double
Definition: G4Types.hh:76
Definition of the B03ImportanceDetectorConstruction class.
const G4VPhysicalVolume * GetPVolume(const G4String &name) const
void DivideByVelocity(G4bool flg=true)
void AddPVolume(const G4GeometryCell &cell)