Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B03ImportanceDetectorConstruction Class Reference

#include <B03ImportanceDetectorConstruction.hh>

Inheritance diagram for B03ImportanceDetectorConstruction:
Collaboration diagram for B03ImportanceDetectorConstruction:

Public Member Functions

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

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 47 of file B03ImportanceDetectorConstruction.hh.

Constructor & Destructor Documentation

B03ImportanceDetectorConstruction::B03ImportanceDetectorConstruction ( G4String  worldName)

Definition at line 58 of file B03ImportanceDetectorConstruction.cc.

59 :G4VUserParallelWorld(worldName),fLogicalVolumeVector()
60 {
61  // Construct();
62 }
G4VUserParallelWorld(G4String worldName)
B03ImportanceDetectorConstruction::~B03ImportanceDetectorConstruction ( )

Definition at line 66 of file B03ImportanceDetectorConstruction.cc.

67 {
68  fLogicalVolumeVector.clear();
69 }

Member Function Documentation

void B03ImportanceDetectorConstruction::ConstructSD ( )
virtual

Reimplemented from G4VUserParallelWorld.

Definition at line 262 of file B03ImportanceDetectorConstruction.cc.

263 {
264 
266  //
267  // Sensitive Detector Name
268  G4String concreteSDname = "ConcreteSD";
269 
270  //------------------------
271  // MultiFunctionalDetector
272  //------------------------
273  //
274  // Define MultiFunctionalDetector with name.
275  G4MultiFunctionalDetector* MFDet =
276  new G4MultiFunctionalDetector(concreteSDname);
277  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
278 
279  G4String fltName,particleName;
280  G4SDParticleFilter* neutronFilter =
281  new G4SDParticleFilter(fltName="neutronFilter", particleName="neutron");
282 
283  MFDet->SetFilter(neutronFilter);
284 
285  for (std::vector<G4LogicalVolume *>::iterator it =
286  fLogicalVolumeVector.begin();
287  it != fLogicalVolumeVector.end(); it++){
288  // (*it)->SetSensitiveDetector(MFDet);
289  SetSensitiveDetector((*it)->GetName(), MFDet);
290  }
291 
292  G4String psName;
293  G4PSNofCollision* scorer0 = new G4PSNofCollision(psName="Collisions");
294  MFDet->RegisterPrimitive(scorer0);
295 
296  G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
297  scorer1->Weighted(true);
298  MFDet->RegisterPrimitive(scorer1);
299 
300  G4PSPopulation* scorer2 = new G4PSPopulation(psName="Population");
301  MFDet->RegisterPrimitive(scorer2);
302 
303  G4PSTrackCounter* scorer3 =
304  new G4PSTrackCounter(psName="TrackEnter",fCurrent_In);
305  MFDet->RegisterPrimitive(scorer3);
306 
307  G4PSTrackLength* scorer4 = new G4PSTrackLength(psName="SL");
308  MFDet->RegisterPrimitive(scorer4);
309 
310  G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
311  scorer5->Weighted(true);
312  MFDet->RegisterPrimitive(scorer5);
313 
314  G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
315  scorer6->Weighted(true);
316  scorer6->MultiplyKineticEnergy(true);
317  MFDet->RegisterPrimitive(scorer6);
318 
319  G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
320  scorer7->Weighted(true);
321  scorer7->DivideByVelocity(true);
322  MFDet->RegisterPrimitive(scorer7);
323 
324  G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
325  scorer8->Weighted(true);
326  scorer8->MultiplyKineticEnergy(true);
327  scorer8->DivideByVelocity(true);
328  MFDet->RegisterPrimitive(scorer8);
329 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
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:

G4String B03ImportanceDetectorConstruction::GetCellName ( G4int  i)

Definition at line 200 of file B03ImportanceDetectorConstruction.cc.

200  {
201  std::ostringstream os;
202  os << "cell_";
203  if (i<10) {
204  os << "0";
205  }
206  os << i;
207  G4String name = os.str();
208  return name;
209 }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

G4GeometryCell B03ImportanceDetectorConstruction::GetGeometryCell ( G4int  i)

Definition at line 213 of file B03ImportanceDetectorConstruction.cc.

213  {
215  const G4VPhysicalVolume *p=0;
216  p = fPVolumeStore.GetPVolume(name);
217  if (p) {
218  return G4GeometryCell(*p,0);
219  }
220  else {
221  G4cout << "B03ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
222  << " couldn't get G4GeometryCell" << G4endl;
223  return G4GeometryCell(*fGhostWorld,-2);
224  }
225 }
const XML_Char * name
Definition: expat.h:151
const char * p
Definition: xmltok.h:285
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
const G4VPhysicalVolume * GetPVolume(const G4String &name) const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 187 of file B03ImportanceDetectorConstruction.cc.

187  {
188  return *fPVolumeStore.GetPVolume(name);
189 }
const G4VPhysicalVolume * GetPVolume(const G4String &name) const

Here is the call graph for this function:

G4VPhysicalVolume * B03ImportanceDetectorConstruction::GetWorldVolume ( )

Definition at line 236 of file B03ImportanceDetectorConstruction.cc.

236  {
237  return fGhostWorld;
238 }

Here is the caller graph for this function:

G4VPhysicalVolume & B03ImportanceDetectorConstruction::GetWorldVolumeAddress ( ) const

Definition at line 230 of file B03ImportanceDetectorConstruction.cc.

230  {
231  return *fGhostWorld;
232 }

Here is the caller graph for this function:

G4String B03ImportanceDetectorConstruction::ListPhysNamesAsG4String ( )

Definition at line 193 of file B03ImportanceDetectorConstruction.cc.

193  {
194  G4String names(fPVolumeStore.GetPNames());
195  return names;
196 }
G4String GetPNames() const

Here is the call graph for this function:

void B03ImportanceDetectorConstruction::SetSensitive ( )

Definition at line 242 of file B03ImportanceDetectorConstruction.cc.

242  {
243 
244  // -------------------------------------------------
245  // The collection names of defined Primitives are
246  // 0 ConcreteSD/Collisions
247  // 1 ConcreteSD/CollWeight
248  // 2 ConcreteSD/Population
249  // 3 ConcreteSD/TrackEnter
250  // 4 ConcreteSD/SL
251  // 5 ConcreteSD/SLW
252  // 6 ConcreteSD/SLWE
253  // 7 ConcreteSD/SLW_V
254  // 8 ConcreteSD/SLWE_V
255  // -------------------------------------------------
256 
257  //now moved to ConstructSD()
258 
259 }

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