Geant4  10.02.p03
B02ImportanceDetectorConstruction Class Reference

#include <B02ImportanceDetectorConstruction.hh>

Inheritance diagram for B02ImportanceDetectorConstruction:
Collaboration diagram for B02ImportanceDetectorConstruction:

Public Member Functions

 B02ImportanceDetectorConstruction (G4String worldName)
 
virtual ~B02ImportanceDetectorConstruction ()
 
const G4VPhysicalVolumeGetPhysicalVolumeByName (const G4String &name) const
 
G4VPhysicalVolumeGetWorldVolumeAddress () const
 
G4String ListPhysNamesAsG4String ()
 
G4String GetCellName (G4int i)
 
G4GeometryCell GetGeometryCell (G4int i)
 
G4VPhysicalVolumeGetWorldVolume ()
 
void SetSensitive ()
 
virtual void Construct ()
 
virtual void ConstructSD ()
 
G4VIStoreCreateImportanceStore ()
 
G4VWeightWindowStoreCreateWeightWindowStore ()
 
- Public Member Functions inherited from G4VUserParallelWorld
 G4VUserParallelWorld (G4String worldName)
 
virtual ~G4VUserParallelWorld ()
 
G4String GetName ()
 

Private Attributes

B02PVolumeStore fPVolumeStore
 
std::vector< G4LogicalVolume *> fLogicalVolumeVector
 
G4VPhysicalVolumefGhostWorld
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserParallelWorld
G4VPhysicalVolumeGetWorld ()
 
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 
- Protected Attributes inherited from G4VUserParallelWorld
G4String fWorldName
 

Detailed Description

Definition at line 51 of file B02ImportanceDetectorConstruction.hh.

Constructor & Destructor Documentation

◆ B02ImportanceDetectorConstruction()

B02ImportanceDetectorConstruction::B02ImportanceDetectorConstruction ( G4String  worldName)

Definition at line 64 of file B02ImportanceDetectorConstruction.cc.

66 {
67  // Construct();
68 }
std::vector< G4LogicalVolume *> fLogicalVolumeVector
G4VUserParallelWorld(G4String worldName)

◆ ~B02ImportanceDetectorConstruction()

B02ImportanceDetectorConstruction::~B02ImportanceDetectorConstruction ( )
virtual

Definition at line 72 of file B02ImportanceDetectorConstruction.cc.

73 {
74  fLogicalVolumeVector.clear();
75 }
std::vector< G4LogicalVolume *> fLogicalVolumeVector

Member Function Documentation

◆ Construct()

void B02ImportanceDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserParallelWorld.

Definition at line 79 of file B02ImportanceDetectorConstruction.cc.

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
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));
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);
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);
189 
190  SetSensitive();
191 
192 }
static const double cm
Definition: G4SIunits.hh:118
CLHEP::Hep3Vector G4ThreeVector
G4VPhysicalVolume * GetWorld()
G4String name
Definition: TRTMaterials.hh:40
Definition: G4Tubs.hh:85
std::vector< G4LogicalVolume *> fLogicalVolumeVector
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static const double deg
Definition: G4SIunits.hh:151
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
void AddPVolume(const G4GeometryCell &cell)
Here is the call graph for this function:

◆ ConstructSD()

void B02ImportanceDetectorConstruction::ConstructSD ( )
virtual

Reimplemented from G4VUserParallelWorld.

Definition at line 272 of file B02ImportanceDetectorConstruction.cc.

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 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
std::vector< G4LogicalVolume *> fLogicalVolumeVector
void Weighted(G4bool flg=true)
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
void Weighted(G4bool flg=true)
void SetFilter(G4VSDFilter *value)
void AddNewDetector(G4VSensitiveDetector *aSD)
Definition: G4SDManager.cc:71
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void MultiplyKineticEnergy(G4bool flg=true)
void DivideByVelocity(G4bool flg=true)
Here is the call graph for this function:

◆ CreateImportanceStore()

G4VIStore * B02ImportanceDetectorConstruction::CreateImportanceStore ( )

Definition at line 346 of file B02ImportanceDetectorConstruction.cc.

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 }
G4int GetReplicaNumber() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void AddImportanceGeometryCell(G4double importance, const G4GeometryCell &gCell)
Definition: G4IStore.cc:103
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
static G4IStore * GetInstance()
Definition: G4IStore.cc:210
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4VPhysicalVolume & GetPhysicalVolume() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateWeightWindowStore()

G4VWeightWindowStore * B02ImportanceDetectorConstruction::CreateWeightWindowStore ( )

Definition at line 403 of file B02ImportanceDetectorConstruction.cc.

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 }
G4int GetReplicaNumber() const
void AddLowerWeights(const G4GeometryCell &gCell, const std::vector< G4double > &lowerWeights)
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
static const double GeV
Definition: G4SIunits.hh:214
static G4WeightWindowStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
void SetGeneralUpperEnergyBounds(const std::set< G4double, std::less< G4double > > &enBounds)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4VPhysicalVolume & GetPhysicalVolume() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCellName()

G4String B02ImportanceDetectorConstruction::GetCellName ( G4int  i)

Definition at line 210 of file B02ImportanceDetectorConstruction.cc.

210  {
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 }
G4String name
Definition: TRTMaterials.hh:40
Here is the caller graph for this function:

◆ GetGeometryCell()

G4GeometryCell B02ImportanceDetectorConstruction::GetGeometryCell ( G4int  i)

Definition at line 223 of file B02ImportanceDetectorConstruction.cc.

223  {
225  const G4VPhysicalVolume *p=0;
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 }
G4String name
Definition: TRTMaterials.hh:40
const G4VPhysicalVolume * GetPVolume(const G4String &name) const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPhysicalVolumeByName()

const G4VPhysicalVolume & B02ImportanceDetectorConstruction::GetPhysicalVolumeByName ( const G4String name) const

Definition at line 197 of file B02ImportanceDetectorConstruction.cc.

197  {
198  return *fPVolumeStore.GetPVolume(name);
199 }
const G4VPhysicalVolume * GetPVolume(const G4String &name) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetWorldVolume()

G4VPhysicalVolume * B02ImportanceDetectorConstruction::GetWorldVolume ( )

Definition at line 246 of file B02ImportanceDetectorConstruction.cc.

246  {
247  return fGhostWorld;
248 }
Here is the caller graph for this function:

◆ GetWorldVolumeAddress()

G4VPhysicalVolume & B02ImportanceDetectorConstruction::GetWorldVolumeAddress ( ) const

Definition at line 240 of file B02ImportanceDetectorConstruction.cc.

240  {
241  return *fGhostWorld;
242 }
Here is the caller graph for this function:

◆ ListPhysNamesAsG4String()

G4String B02ImportanceDetectorConstruction::ListPhysNamesAsG4String ( )

Definition at line 203 of file B02ImportanceDetectorConstruction.cc.

Here is the call graph for this function:

◆ SetSensitive()

void B02ImportanceDetectorConstruction::SetSensitive ( )

Definition at line 252 of file B02ImportanceDetectorConstruction.cc.

252  {
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 }
Here is the caller graph for this function:

Member Data Documentation

◆ fGhostWorld

G4VPhysicalVolume* B02ImportanceDetectorConstruction::fGhostWorld
private

Definition at line 88 of file B02ImportanceDetectorConstruction.hh.

◆ fLogicalVolumeVector

std::vector< G4LogicalVolume * > B02ImportanceDetectorConstruction::fLogicalVolumeVector
private

Definition at line 84 of file B02ImportanceDetectorConstruction.hh.

◆ fPVolumeStore

B02PVolumeStore B02ImportanceDetectorConstruction::fPVolumeStore
private

Definition at line 79 of file B02ImportanceDetectorConstruction.hh.


The documentation for this class was generated from the following files: