Geant4  10.00.p02
RE02DetectorConstruction Class Reference

Uer detector construction class. More...

#include <RE02DetectorConstruction.hh>

+ Inheritance diagram for RE02DetectorConstruction:
+ Collaboration diagram for RE02DetectorConstruction:

Public Member Functions

 RE02DetectorConstruction ()
 
virtual ~RE02DetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
virtual void ConstructSDandField ()
 
void SetPhantomSize (G4ThreeVector size)
 
const G4ThreeVectorGetPhantomSize () const
 
void SetNumberOfSegmentsInPhantom (G4int nx, G4int ny, G4int nz)
 
void GetNumberOfSegmentsInPhantom (G4int &nx, G4int &ny, G4int &nz) const
 
void SetLeadSegment (G4bool flag=TRUE)
 
G4bool IsLeadSegment ()
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Private Attributes

G4ThreeVector fPhantomSize
 
G4int fNx
 
G4int fNy
 
G4int fNz
 
G4bool fInsertLead
 
G4LogicalVolumefLVPhantomSens
 

Additional Inherited Members

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

Detailed Description

Uer detector construction class.

(Description)

Detector construction for example RE02.

[Geometry] The world volume is defined as 200 cm x 200 cm x 200 cm box with Air. Water phantom is defined as 200 mm x 200 mm x 400 mm box with Water. The water phantom is divided into 100 segments in x,y plane using replication, and then divided into 200 segments perpendicular to z axis using nested parameterised volume. These values are defined at constructor, e.g. the size of water phantom (fPhantomSize), and number of segmentation of water phantom (fNx, fNy, fNz).

By default, lead plates are inserted into the position of even order segments. NIST database is used for materials.

[Scorer] Assignment of G4MultiFunctionalDetector and G4PrimitiveScorer

is demonstrated in this example.

The collection names of defined Primitives are 0 PhantomSD/totalEDep 1 PhantomSD/protonEDep 2 PhantomSD/protonNStep 3 PhantomSD/chargedPassCellFlux 4 PhantomSD/chargedCellFlux 5 PhantomSD/chargedSurfFlux 6 PhantomSD/gammaSurfCurr000 7 PhantomSD/gammaSurfCurr001 9 PhantomSD/gammaSurdCurr002

10 PhantomSD/gammaSurdCurr003

Please see README for detail description.

Definition at line 111 of file RE02DetectorConstruction.hh.

Constructor & Destructor Documentation

RE02DetectorConstruction::RE02DetectorConstruction ( )

Definition at line 106 of file RE02DetectorConstruction.cc.

References fInsertLead, fNx, fNy, fNz, fPhantomSize, mm, and TRUE.

RE02DetectorConstruction::~RE02DetectorConstruction ( )
virtual

Definition at line 118 of file RE02DetectorConstruction.cc.

Member Function Documentation

G4VPhysicalVolume * RE02DetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 122 of file RE02DetectorConstruction.cc.

References cm, G4NistManager::FindOrBuildMaterial(), fLVPhantomSens, fNx, fNy, fNz, fPhantomSize, G4cout, G4endl, G4Material::GetMaterialTable(), G4NistManager::Instance(), G4VisAttributes::Invisible, IsLeadSegment(), kUndefined, kXAxis, kYAxis, mm, and G4LogicalVolume::SetVisAttributes().

+ Here is the call graph for this function:

void RE02DetectorConstruction::ConstructSDandField ( )
virtual
void RE02DetectorConstruction::GetNumberOfSegmentsInPhantom ( G4int nx,
G4int ny,
G4int nz 
) const
inline

Definition at line 131 of file RE02DetectorConstruction.hh.

References fNx, fNy, and fNz.

Referenced by RE02RunAction::EndOfRunAction().

+ Here is the caller graph for this function:

const G4ThreeVector& RE02DetectorConstruction::GetPhantomSize ( ) const
inline

Definition at line 127 of file RE02DetectorConstruction.hh.

References fPhantomSize.

G4bool RE02DetectorConstruction::IsLeadSegment ( )
inline

Definition at line 135 of file RE02DetectorConstruction.hh.

References fInsertLead.

Referenced by Construct().

+ Here is the caller graph for this function:

void RE02DetectorConstruction::SetLeadSegment ( G4bool  flag = TRUE)
inline

Definition at line 134 of file RE02DetectorConstruction.hh.

References fInsertLead.

Referenced by main().

+ Here is the caller graph for this function:

void RE02DetectorConstruction::SetNumberOfSegmentsInPhantom ( G4int  nx,
G4int  ny,
G4int  nz 
)
inline

Definition at line 129 of file RE02DetectorConstruction.hh.

References fNx, fNy, and fNz.

Referenced by main().

+ Here is the caller graph for this function:

void RE02DetectorConstruction::SetPhantomSize ( G4ThreeVector  size)
inline

Definition at line 126 of file RE02DetectorConstruction.hh.

References fPhantomSize.

Referenced by main().

+ Here is the caller graph for this function:

Member Data Documentation

G4bool RE02DetectorConstruction::fInsertLead
private
G4LogicalVolume* RE02DetectorConstruction::fLVPhantomSens
private

Definition at line 142 of file RE02DetectorConstruction.hh.

Referenced by Construct(), and ConstructSDandField().

G4ThreeVector RE02DetectorConstruction::fPhantomSize
private

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