Geant4  10.02.p03
DicomRegularDetectorConstruction Class Reference

#include <DicomRegularDetectorConstruction.hh>

Inheritance diagram for DicomRegularDetectorConstruction:
Collaboration diagram for DicomRegularDetectorConstruction:

Public Member Functions

 DicomRegularDetectorConstruction ()
 
 ~DicomRegularDetectorConstruction ()
 
- Public Member Functions inherited from DicomDetectorConstruction
 DicomDetectorConstruction ()
 
 ~DicomDetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
- 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 Member Functions

virtual void ConstructPhantom ()
 

Additional Inherited Members

- Protected Member Functions inherited from DicomDetectorConstruction
void InitialisationOfMaterials ()
 
void ReadPhantomData ()
 
void ReadPhantomDataFile (const G4String &fname)
 
void MergeZSliceHeaders ()
 
G4MaterialBuildMaterialWithChangingDensity (const G4Material *origMate, float density, G4String newMateName)
 
void ConstructPhantomContainer ()
 
void SetScorer (G4LogicalVolume *voxel_logic)
 
virtual void ConstructSDandField ()
 
- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 
- Protected Attributes inherited from DicomDetectorConstruction
G4MaterialfAir
 
G4BoxfWorld_solid
 
G4LogicalVolumefWorld_logic
 
G4VPhysicalVolumefWorld_phys
 
G4BoxfContainer_solid
 
G4LogicalVolumefContainer_logic
 
G4VPhysicalVolumefContainer_phys
 
G4int fNoFiles
 
std::vector< G4Material * > fOriginalMaterials
 
std::vector< G4Material * > fMaterials
 
size_t * fMateIDs
 
std::map< G4int, G4doublefDensityDiffs
 
std::vector< DicomPhantomZSliceHeader * > fZSliceHeaders
 
DicomPhantomZSliceHeaderfZSliceHeaderMerged
 
G4int fNVoxelX
 
G4int fNVoxelY
 
G4int fNVoxelZ
 
G4double fVoxelHalfDimX
 
G4double fVoxelHalfDimY
 
G4double fVoxelHalfDimZ
 
DicomPhantomZSliceMergedfMergedSlices
 
std::set< G4LogicalVolume * > fScorers
 
G4bool fConstructed
 

Detailed Description

DicomRegularDetectorConstruction class

Construct the phantom using DicomPhantomParameterisatin

History: 30.11.07 First version

Author
P. Arce

Definition at line 51 of file DicomRegularDetectorConstruction.hh.

Constructor & Destructor Documentation

◆ DicomRegularDetectorConstruction()

DicomRegularDetectorConstruction::DicomRegularDetectorConstruction ( )

◆ ~DicomRegularDetectorConstruction()

DicomRegularDetectorConstruction::~DicomRegularDetectorConstruction ( )

Definition at line 59 of file DicomRegularDetectorConstruction.cc.

60 {
61 }

Member Function Documentation

◆ ConstructPhantom()

void DicomRegularDetectorConstruction::ConstructPhantom ( )
privatevirtual

Implements DicomDetectorConstruction.

Definition at line 64 of file DicomRegularDetectorConstruction.cc.

65 {
66 #ifdef G4VERBOSE
67  G4cout << "DicomRegularDetectorConstruction::ConstructPhantom " << G4endl;
68 #endif
69 
70  //----- Create parameterisation
73 
74  //----- Set voxel dimensions
76 
77  //----- Set number of voxels
79 
80  //----- Set list of materials
81  param->SetMaterials( fMaterials );
82 
83  //----- Set list of material indices: for each voxel it is a number that
84  // correspond to the index of its material in the vector of materials
85  // defined above
86  param->SetMaterialIndices( fMateIDs );
87 
88  //----- Define voxel logical volume
89  G4Box* voxel_solid =
91  G4LogicalVolume* voxel_logic =
92  new G4LogicalVolume(voxel_solid,fMaterials[0],"VoxelLogical",
93  0,0,0);
94  // material is not relevant, it will be changed by the
95  // ComputeMaterial method of the parameterisation
96 
97  voxel_logic->SetVisAttributes(
99 
100  //--- Assign the fContainer volume of the parameterisation
102 
103  //--- Assure yourself that the voxels are completely filling the
104  // fContainer volume
108 
109  //----- The G4PVParameterised object that uses the created parameterisation
110  // should be placed in the fContainer logical volume
111  G4PVParameterised * phantom_phys =
112  new G4PVParameterised("phantom",voxel_logic,fContainer_logic,
113  kXAxis, fNVoxelX*fNVoxelY*fNVoxelZ, param);
114  // if axis is set as kUndefined instead of kXAxis, GEANT4 will
115  // do an smart voxel optimisation
116  // (not needed if G4RegularNavigation is used)
117 
118  //----- Set this physical volume as having a regular structure of type 1,
119  // so that G4RegularNavigation is used
120  phantom_phys->SetRegularStructureId(1); // if not set, G4VoxelNavigation
121  //will be used instead
122 
123  SetScorer(voxel_logic);
124 }
void CheckVoxelsFillContainer(G4double contX, G4double contY, G4double contZ) const
void SetScorer(G4LogicalVolume *voxel_logic)
void SetMaterials(std::vector< G4Material *> &mates)
Definition: G4Box.hh:64
void BuildContainerSolid(G4VPhysicalVolume *pPhysicalVol)
G4double GetXHalfLength() const
G4GLOB_DLL std::ostream G4cout
void SetVoxelDimensions(G4double halfx, G4double halfy, G4double halfz)
G4double GetZHalfLength() const
G4double GetYHalfLength() const
virtual void SetRegularStructureId(G4int Code)
void SetNoVoxel(size_t nx, size_t ny, size_t nz)
static const G4VisAttributes Invisible
#define G4endl
Definition: G4ios.hh:61
void SetMaterialIndices(size_t *matInd)
std::vector< G4Material * > fMaterials
Class inherited from G4PhantomParameterisation to provide different.
void SetVisAttributes(const G4VisAttributes *pVA)
Here is the call graph for this function:

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