Geant4  10.03.p03
 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: B02ImportanceDetectorConstruction.cc 77336 2013-11-22 15:04:18Z ahoward$
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 // for importance biasing
56 #include "G4IStore.hh"
57 
58 // for weight window technique
59 #include "G4WeightWindowStore.hh"
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
65 :G4VUserParallelWorld(worldName),fLogicalVolumeVector()
66 {
67  // Construct();
68 }
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 {
74  fLogicalVolumeVector.clear();
75 }
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78 
80 {
81  G4cout << " constructing parallel world " << G4endl;
82 
83  G4Material* dummyMat = 0;
84 
85  //GetWorld methods create a clone of the mass world to the parallel world (!)
86  // via the transportation manager
87  fGhostWorld = GetWorld();
88  G4cout << " B02ImportanceDetectorConstruction:: ghostWorldName = "
89  << fGhostWorld->GetName() << G4endl;
90  G4LogicalVolume* worldLogical = fGhostWorld->GetLogicalVolume();
91  fLogicalVolumeVector.push_back(worldLogical);
92 
93  // fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
94  fPVolumeStore.AddPVolume(G4GeometryCell(*fGhostWorld, 0));
95 
96  // creating 18 slobs of 10 cm thicknes
97 
98  G4double innerRadiusShield = 0*cm;
99  G4double outerRadiusShield = 100*cm;
100  G4double heightShield = 5*cm;
101  G4double startAngleShield = 0*deg;
102  G4double spanningAngleShield = 360*deg;
103 
104  G4Tubs *aShield = new G4Tubs("aShield",
105  innerRadiusShield,
106  outerRadiusShield,
107  heightShield,
108  startAngleShield,
109  spanningAngleShield);
110 
111  // logical parallel cells
112 
113  G4LogicalVolume *aShield_log_imp =
114  new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
115  fLogicalVolumeVector.push_back(aShield_log_imp);
116 
117  // physical parallel cells
118  G4String name = "none";
119  G4int i = 1;
120  G4double startz = -85*cm;
121  // for (i=1; i<=18; ++i) {
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*heightShield);
129  G4VPhysicalVolume *pvol =
130  new G4PVPlacement(0,
131  G4ThreeVector(pos_x, pos_y, pos_z),
132  aShield_log_imp,
133  name,
134  worldLogical,
135  false,
136  i);
137  // 0);
138  G4GeometryCell cell(*pvol, i);
139  // G4GeometryCell cell(*pvol, 0);
140  fPVolumeStore.AddPVolume(cell);
141  }
142 
143  // filling the rest of the world volumr behind the concrete with
144  // another slob which should get the same importance value as the
145  // last slob
146  innerRadiusShield = 0*cm;
147  // outerRadiusShield = 110*cm; exceeds world volume!!!!
148  outerRadiusShield = 100*cm;
149  // heightShield = 10*cm;
150  heightShield = 5*cm;
151  startAngleShield = 0*deg;
152  spanningAngleShield = 360*deg;
153 
154  G4Tubs *aRest = new G4Tubs("Rest",
155  innerRadiusShield,
156  outerRadiusShield,
157  heightShield,
158  startAngleShield,
159  spanningAngleShield);
160 
161  G4LogicalVolume *aRest_log =
162  new G4LogicalVolume(aRest, dummyMat, "aRest_log");
163 
164  fLogicalVolumeVector.push_back(aRest_log);
165 
166  name = GetCellName(19);
167 
168  G4double pos_x = 0*cm;
169  G4double pos_y = 0*cm;
170  // G4double pos_z = 100*cm;
171  G4double pos_z = 95*cm;
172  G4VPhysicalVolume *pvol =
173  new G4PVPlacement(0,
174  G4ThreeVector(pos_x, pos_y, pos_z),
175  aRest_log,
176  name,
177  worldLogical,
178  false,
179  19);
180  // 0);
181  G4GeometryCell cell(*pvol, 19);
182  // G4GeometryCell cell(*pvol, 0);
183  fPVolumeStore.AddPVolume(cell);
184 
185  SetSensitive();
186 
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190 
193  return *fPVolumeStore.GetPVolume(name);
194 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197 
199  G4String names(fPVolumeStore.GetPNames());
200  return names;
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204 
206  std::ostringstream os;
207  os << "cell_";
208  if (i<10) {
209  os << "0";
210  }
211  os << i;
212  G4String name = os.str();
213  return name;
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
220  const G4VPhysicalVolume *p=0;
221  p = fPVolumeStore.GetPVolume(name);
222  if (p) {
223  return G4GeometryCell(*p,0);
224  }
225  else {
226  G4cout << "B02ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
227  << " couldn't get G4GeometryCell" << G4endl;
228  return G4GeometryCell(*fGhostWorld,-2);
229  }
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
233 
236  return *fGhostWorld;
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240 
242  return fGhostWorld;
243 }
244 
245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
246 
248 
249  // -------------------------------------------------
250  // The collection names of defined Primitives are
251  // 0 ConcreteSD/Collisions
252  // 1 ConcreteSD/CollWeight
253  // 2 ConcreteSD/Population
254  // 3 ConcreteSD/TrackEnter
255  // 4 ConcreteSD/SL
256  // 5 ConcreteSD/SLW
257  // 6 ConcreteSD/SLWE
258  // 7 ConcreteSD/SLW_V
259  // 8 ConcreteSD/SLWE_V
260  // -------------------------------------------------
261 
262  // moved to ConstructSD() for MT compliance
263 
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
268 {
269 
271  //
272  // Sensitive Detector Name
273  G4String concreteSDname = "ConcreteSD";
274 
275  //------------------------
276  // MultiFunctionalDetector
277  //------------------------
278  //
279  // Define MultiFunctionalDetector with name.
280  G4MultiFunctionalDetector* MFDet =
281  new G4MultiFunctionalDetector(concreteSDname);
282  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
283 
284  G4String fltName,particleName;
285  G4SDParticleFilter* neutronFilter =
286  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
287 
288  MFDet->SetFilter(neutronFilter);
289 
290  for (std::vector<G4LogicalVolume *>::iterator it =
291  fLogicalVolumeVector.begin();
292  it != fLogicalVolumeVector.end(); it++){
293  // (*it)->SetSensitiveDetector(MFDet);
294  SetSensitiveDetector((*it)->GetName(), MFDet);
295  }
296 
297  G4String psName;
298  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
299  MFDet->RegisterPrimitive(scorer0);
300 
301  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
302  scorer1->Weighted(true);
303  MFDet->RegisterPrimitive(scorer1);
304 
305  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
306  MFDet->RegisterPrimitive(scorer2);
307 
308  G4PSTrackCounter* scorer3 =
309  new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
310  MFDet->RegisterPrimitive(scorer3);
311 
312  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
313  MFDet->RegisterPrimitive(scorer4);
314 
315  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
316  scorer5->Weighted(true);
317  MFDet->RegisterPrimitive(scorer5);
318 
319  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
320  scorer6->Weighted(true);
321  scorer6->MultiplyKineticEnergy(true);
322  MFDet->RegisterPrimitive(scorer6);
323 
324  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
325  scorer7->Weighted(true);
326  scorer7->DivideByVelocity(true);
327  MFDet->RegisterPrimitive(scorer7);
328 
329  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
330  scorer8->Weighted(true);
331  scorer8->MultiplyKineticEnergy(true);
332  scorer8->DivideByVelocity(true);
333  MFDet->RegisterPrimitive(scorer8);
334 }
335 
336 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
338 {
339 
340  G4cout << " B02ImportanceDetectorConstruction:: Creating Importance Store "
341  << G4endl;
342  if (!fPVolumeStore.Size())
343  {
344  G4Exception("B02ImportanceDetectorConstruction::CreateImportanceStore"
345  ,"exampleB02_0001",RunMustBeAborted
346  ,"no physical volumes created yet!");
347  }
348 
349  // creating and filling the importance store
350 
351  // G4IStore *istore = new G4IStore(*fWorldVolume);
352 
354 
355  G4GeometryCell gWorldVolumeCell(GetWorldVolumeAddress(), 0);
356 
357  G4double imp =1;
358 
359  istore->AddImportanceGeometryCell(1, gWorldVolumeCell);
360 
361  // set importance values and create scorers
362  G4int cell(1);
363  for (cell=1; cell<=18; cell++) {
364  G4GeometryCell gCell = GetGeometryCell(cell);
365  G4cout << " adding cell: " << cell
366  << " replica: " << gCell.GetReplicaNumber()
367  << " name: " << gCell.GetPhysicalVolume().GetName() << G4endl;
368  imp = std::pow(2.0,cell-1);
369 
370  G4cout << "Going to assign importance: " << imp << ", to volume: "
371  << gCell.GetPhysicalVolume().GetName() << G4endl;
372  //x aIstore.AddImportanceGeometryCell(imp, gCell);
373  istore->AddImportanceGeometryCell(imp, gCell.GetPhysicalVolume(), cell);
374  }
375 
376  // creating the geometry cell and add both to the store
377  // G4GeometryCell gCell = GetGeometryCell(18);
378 
379  // create importance geometry cell pair for the "rest"cell
380  // with the same importance as the last concrete cell
381  G4GeometryCell gCell = GetGeometryCell(19);
382  // G4double imp = std::pow(2.0,18);
383  imp = std::pow(2.0,17);
384  istore->AddImportanceGeometryCell(imp, gCell.GetPhysicalVolume(), 19);
385 
386  return istore;
387 
388 }
389 
390 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
391 
394 {
395  G4cout << " B02ImportanceDetectorConstruction:: Creating Importance Store "
396  << G4endl;
397  if (!fPVolumeStore.Size())
398  {
399  G4Exception("B02ImportanceDetectorConstruction::CreateWeightWindowStore"
400  ,"exampleB02_0002",RunMustBeAborted
401  ,"no physical volumes created yet!");
402  }
403 
404  // creating and filling the importance store
405 
406  // G4IStore *istore = new G4IStore(*fWorldVolume);
407 
409 
410  // create one energy region covering the energies of the problem
411  //
412  std::set<G4double, std::less<G4double> > enBounds;
413  enBounds.insert(1 * GeV);
414  wwstore->SetGeneralUpperEnergyBounds(enBounds);
415 
416  G4int n = 0;
417  G4double lowerWeight =1;
418  std::vector<G4double> lowerWeights;
419 
420  lowerWeights.push_back(1);
421  G4GeometryCell gWorldCell(GetWorldVolumeAddress(),0);
422  wwstore->AddLowerWeights(gWorldCell, lowerWeights);
423 
424  G4int cell(1);
425  for (cell=1; cell<=18; cell++) {
426  G4GeometryCell gCell = GetGeometryCell(cell);
427  G4cout << " adding cell: " << cell
428  << " replica: " << gCell.GetReplicaNumber()
429  << " name: " << gCell.GetPhysicalVolume().GetName() << G4endl;
430 
431  lowerWeight = 1./std::pow(2., n++);
432  G4cout << "Going to assign lower weight: " << lowerWeight
433  << ", to volume: "
434  << gCell.GetPhysicalVolume().GetName() << G4endl;
435  lowerWeights.clear();
436  lowerWeights.push_back(lowerWeight);
437  wwstore->AddLowerWeights(gCell, lowerWeights);
438  }
439 
440  // the remaining part pf the geometry (rest) gets the same
441  // lower weight bound as the last conrete cell
442  //
443 
444  // create importance geometry cell pair for the "rest"cell
445  // with the same importance as the last concrete cell
446  G4GeometryCell gCell = GetGeometryCell(19);
447  wwstore->AddLowerWeights(gCell, lowerWeights);
448 
449  return wwstore;
450 
451 }
const XML_Char * name
Definition: expat.h:151
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
CLHEP::Hep3Vector G4ThreeVector
G4VPhysicalVolume * GetWorld()
const char * p
Definition: xmltok.h:285
void AddLowerWeights(const G4GeometryCell &gCell, const std::vector< G4double > &lowerWeights)
Definition: G4Tubs.hh:85
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
static constexpr double cm
Definition: G4SIunits.hh:119
G4int GetReplicaNumber() const
const G4VPhysicalVolume * GetPVolume(const G4String &name) const
void AddImportanceGeometryCell(G4double importance, const G4GeometryCell &gCell)
Definition: G4IStore.cc:107
void SetFilter(G4VSDFilter *value)
const G4int n
static G4WeightWindowStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
const G4VPhysicalVolume & GetPhysicalVolume() const
G4String GetPNames() const
static G4IStore * GetInstance()
Definition: G4IStore.cc:236
G4LogicalVolume * GetLogicalVolume() const
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void MultiplyKineticEnergy(G4bool flg=true)
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
static constexpr double deg
Definition: G4SIunits.hh:152
const G4VPhysicalVolume & GetPhysicalVolumeByName(const G4String &name) const
void DivideByVelocity(G4bool flg=true)
Definition of the B02ImportanceDetectorConstruction class.
void AddPVolume(const G4GeometryCell &cell)