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