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

#include <HadrontherapyDetectorConstruction.hh>

Collaboration diagram for HadrontherapyDetectorConstruction:

Public Member Functions

 HadrontherapyDetectorConstruction (G4VPhysicalVolume *)
 
 ~HadrontherapyDetectorConstruction ()
 
void InitializeDetectorROGeometry (HadrontherapyDetectorROGeometry *, G4ThreeVector detectorToWorldPosition)
 
G4ThreeVector GetDetectorToWorldPosition ()
 
G4ThreeVector GetDetectorToPhantomPosition ()
 
void SetDetectorPosition ()
 
bool IsInside (G4double detectorX, G4double detectorY, G4double detectorZ, G4double phantomX, G4double phantomY, G4double phantomZ, G4ThreeVector pos)
 
G4bool SetPhantomMaterial (G4String material)
 
void SetVoxelSize (G4double sizeX, G4double sizeY, G4double sizeZ)
 
void SetDetectorSize (G4double sizeX, G4double sizeY, G4double sizeZ)
 
void SetPhantomSize (G4double sizeX, G4double sizeY, G4double sizeZ)
 
void SetPhantomPosition (G4ThreeVector)
 
void SetDetectorToPhantomPosition (G4ThreeVector DetectorToPhantomPosition)
 
void UpdateGeometry ()
 
void PrintParameters ()
 
G4LogicalVolumeGetDetectorLogicalVolume ()
 

Static Public Member Functions

static
HadrontherapyDetectorConstruction
GetInstance ()
 

Public Attributes

G4VPhysicalVolumemotherPhys
 
HadrontherapyDetectorSDdetectorSD
 

Detailed Description

Definition at line 48 of file HadrontherapyDetectorConstruction.hh.

Constructor & Destructor Documentation

HadrontherapyDetectorConstruction::HadrontherapyDetectorConstruction ( G4VPhysicalVolume physicalTreatmentRoom)

Definition at line 62 of file HadrontherapyDetectorConstruction.cc.

63 : motherPhys(physicalTreatmentRoom), // pointer to WORLD volume
64 detectorSD(0), detectorROGeometry(0), matrix(0),
65 phantom(0), detector(0),
66 phantomLogicalVolume(0), detectorLogicalVolume(0),
67 phantomPhysicalVolume(0), detectorPhysicalVolume(0),
68 aRegion(0)
69 {
71 
72  /* NOTE! that the HadrontherapyDetectorConstruction class
73  * does NOT inherit from G4VUserDetectorConstruction G4 class
74  * So the Construct() mandatory virtual method is inside another geometric class
75  * like the passiveProtonBeamLIne, ...
76  */
77 
78  // Messenger to change parameters of the phantom/detector geometry
79  detectorMessenger = new HadrontherapyDetectorMessenger(this);
80 
81  // Default detector voxels size
82  // 200 slabs along the beam direction (X)
83  sizeOfVoxelAlongX = 200 *um;
84  sizeOfVoxelAlongY = 4 *cm;
85  sizeOfVoxelAlongZ = 4 *cm;
86 
87  // Define here the material of the water phantom and of the detector
88  SetPhantomMaterial("G4_WATER");
89  // Construct geometry (messenger commands)
90  SetDetectorSize(4.*cm, 4.*cm, 4.*cm);
91  SetPhantomSize(40. *cm, 40. *cm, 40. *cm);
92  SetPhantomPosition(G4ThreeVector(20. *cm, 0. *cm, 0. *cm));
95  //GetDetectorToWorldPosition();
96 
97  // Write virtual parameters to the real ones and check for consistency
99 }
static HadrontherapyAnalysisManager * GetInstance()
void SetPhantomSize(G4double sizeX, G4double sizeY, G4double sizeZ)
CLHEP::Hep3Vector G4ThreeVector
void SetDetectorSize(G4double sizeX, G4double sizeY, G4double sizeZ)
static constexpr double um
Definition: G4SIunits.hh:113
static constexpr double cm
Definition: G4SIunits.hh:119
void SetDetectorToPhantomPosition(G4ThreeVector DetectorToPhantomPosition)

Here is the call graph for this function:

HadrontherapyDetectorConstruction::~HadrontherapyDetectorConstruction ( )

Definition at line 102 of file HadrontherapyDetectorConstruction.cc.

103 {
104  delete detectorROGeometry;
105  delete matrix;
106  delete detectorMessenger;
107 }

Member Function Documentation

G4LogicalVolume* HadrontherapyDetectorConstruction::GetDetectorLogicalVolume ( )
inline

Definition at line 160 of file HadrontherapyDetectorConstruction.hh.

160 { return detectorLogicalVolume;}
G4ThreeVector HadrontherapyDetectorConstruction::GetDetectorToPhantomPosition ( )
inline

Definition at line 78 of file HadrontherapyDetectorConstruction.hh.

79 {
80  return G4ThreeVector(phantomSizeX/2 - detectorSizeX/2 + detectorPosition.getX(),
81  phantomSizeY/2 - detectorSizeY/2 + detectorPosition.getY(),
82  phantomSizeZ/2 - detectorSizeZ/2 + detectorPosition.getZ()
83  );
84 }
CLHEP::Hep3Vector G4ThreeVector
double getY() const
double getX() const
double getZ() const

Here is the call graph for this function:

G4ThreeVector HadrontherapyDetectorConstruction::GetDetectorToWorldPosition ( )
inline

Definition at line 72 of file HadrontherapyDetectorConstruction.hh.

73  {
74  return phantomPosition + detectorPosition;
75  }

Here is the caller graph for this function:

HadrontherapyDetectorConstruction * HadrontherapyDetectorConstruction::GetInstance ( )
static

Definition at line 110 of file HadrontherapyDetectorConstruction.cc.

111 {
112  return instance;
113 }
static MCTruthManager * instance
void HadrontherapyDetectorConstruction::InitializeDetectorROGeometry ( HadrontherapyDetectorROGeometry RO,
G4ThreeVector  detectorToWorldPosition 
)

Definition at line 213 of file HadrontherapyDetectorConstruction.cc.

217 {
218  RO->Initialize(detectorToWorldPosition,
219  detectorSizeX/2,
220  detectorSizeY/2,
221  detectorSizeZ/2,
222  numberOfVoxelsAlongX,
223  numberOfVoxelsAlongY,
224  numberOfVoxelsAlongZ);
225 }
void Initialize(G4ThreeVector detectorPos, G4double detectorDimX, G4double detectorDimY, G4double detectorDimZ, G4int numberOfVoxelsX, G4int numberOfVoxelsY, G4int numberOfVoxelsZ)

Here is the call graph for this function:

Here is the caller graph for this function:

bool HadrontherapyDetectorConstruction::IsInside ( G4double  detectorX,
G4double  detectorY,
G4double  detectorZ,
G4double  phantomX,
G4double  phantomY,
G4double  phantomZ,
G4ThreeVector  pos 
)
inline

Definition at line 99 of file HadrontherapyDetectorConstruction.hh.

106 {
107 // Dimensions check... X Y and Z
108 // Firstly check what dimension we are modifying
109  {
110  if (detectorX > phantomX)
111  {
112  G4cout << "Error: Detector X dimension must be smaller or equal to the corrispondent of the phantom" << G4endl;
113  return false;
114  }
115  if ( (phantomX - detectorX) < pos.getX())
116  {
117  G4cout << "Error: X dimension doesn't fit with detector to phantom relative position" << G4endl;
118  return false;
119  }
120  }
121 
122  {
123  if (detectorY > phantomY)
124  {
125  G4cout << "Error: Detector Y dimension must be smaller or equal to the corrispondent of the phantom" << G4endl;
126  return false;
127  }
128  if ( (phantomY - detectorY) < pos.getY())
129  {
130  G4cout << "Error: Y dimension doesn't fit with detector to phantom relative position" << G4endl;
131  return false;
132  }
133  }
134 
135  {
136  if (detectorZ > phantomZ)
137  {
138  G4cout << "Error: Detector Z dimension must be smaller or equal to the corrispondent of the phantom" << G4endl;
139  return false;
140  }
141  if ( (phantomZ - detectorZ) < pos.getZ())
142  {
143  G4cout << "Error: Z dimension doesn't fit with detector to phantom relative position" << G4endl;
144  return false;
145  }
146  }
147 
148  return true;
149 }
double getY() const
double getX() const
G4GLOB_DLL std::ostream G4cout
double getZ() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void HadrontherapyDetectorConstruction::PrintParameters ( )

Definition at line 426 of file HadrontherapyDetectorConstruction.cc.

427 {
428 
429  G4cout << "The (X,Y,Z) dimensions of the phantom are : (" <<
430  G4BestUnit( phantom -> GetXHalfLength()*2., "Length") << ',' <<
431  G4BestUnit( phantom -> GetYHalfLength()*2., "Length") << ',' <<
432  G4BestUnit( phantom -> GetZHalfLength()*2., "Length") << ')' << G4endl;
433 
434  G4cout << "The (X,Y,Z) dimensions of the detector are : (" <<
435  G4BestUnit( detector -> GetXHalfLength()*2., "Length") << ',' <<
436  G4BestUnit( detector -> GetYHalfLength()*2., "Length") << ',' <<
437  G4BestUnit( detector -> GetZHalfLength()*2., "Length") << ')' << G4endl;
438 
439  G4cout << "Displacement between Phantom and World is: ";
440  G4cout << "DX= "<< G4BestUnit(phantomPosition.getX(),"Length") <<
441  "DY= "<< G4BestUnit(phantomPosition.getY(),"Length") <<
442  "DZ= "<< G4BestUnit(phantomPosition.getZ(),"Length") << G4endl;
443 
444  G4cout << "The (X,Y,Z) sizes of the Voxels are: (" <<
445  G4BestUnit(sizeOfVoxelAlongX, "Length") << ',' <<
446  G4BestUnit(sizeOfVoxelAlongY, "Length") << ',' <<
447  G4BestUnit(sizeOfVoxelAlongZ, "Length") << ')' << G4endl;
448 
449  G4cout << "The number of Voxels along (X,Y,Z) is: (" <<
450  numberOfVoxelsAlongX << ',' <<
451  numberOfVoxelsAlongY <<',' <<
452  numberOfVoxelsAlongZ << ')' << G4endl;
453 }
double getY() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
double getX() const
G4GLOB_DLL std::ostream G4cout
double getZ() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetDetectorPosition ( )
inline

Definition at line 88 of file HadrontherapyDetectorConstruction.hh.

89  {
90  // Adjust detector position
91  detectorPosition.setX(detectorToPhantomPosition.getX() - phantomSizeX/2 + detectorSizeX/2);
92  detectorPosition.setY(detectorToPhantomPosition.getY() - phantomSizeY/2 + detectorSizeY/2);
93  detectorPosition.setZ(detectorToPhantomPosition.getZ() - phantomSizeZ/2 + detectorSizeZ/2);
94 
95 
96  }
double getY() const
void setY(double)
void setZ(double)
void setX(double)
double getX() const
double getZ() const

Here is the call graph for this function:

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetDetectorSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 291 of file HadrontherapyDetectorConstruction.cc.

292 {
293  if (sizeX > 0.) {detectorSizeX = sizeX;}
294  if (sizeY > 0.) {detectorSizeY = sizeY;}
295  if (sizeZ > 0.) {detectorSizeZ = sizeZ;}
296  SetVoxelSize(sizeOfVoxelAlongX, sizeOfVoxelAlongY, sizeOfVoxelAlongZ);
297 }
void SetVoxelSize(G4double sizeX, G4double sizeY, G4double sizeZ)

Here is the call graph for this function:

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetDetectorToPhantomPosition ( G4ThreeVector  DetectorToPhantomPosition)

Definition at line 314 of file HadrontherapyDetectorConstruction.cc.

315 {
316  detectorToPhantomPosition = displ;
317 }

Here is the caller graph for this function:

G4bool HadrontherapyDetectorConstruction::SetPhantomMaterial ( G4String  material)

Definition at line 255 of file HadrontherapyDetectorConstruction.cc.

256 {
257 
258  if (G4Material* pMat = G4NistManager::Instance()->FindOrBuildMaterial(material, false) )
259  {
260  phantomMaterial = pMat;
261  detectorMaterial = pMat;
262  if (detectorLogicalVolume && phantomLogicalVolume)
263  {
264  detectorLogicalVolume -> SetMaterial(pMat);
265  phantomLogicalVolume -> SetMaterial(pMat);
266 
267  G4RunManager::GetRunManager() -> PhysicsHasBeenModified();
269  G4cout << "The material of Phantom/Detector has been changed to " << material << G4endl;
270  }
271  }
272  else
273  {
274  G4cout << "WARNING: material \"" << material << "\" doesn't exist in NIST elements/materials"
275  " table [located in $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl;
276  G4cout << "Use command \"/parameter/nist\" to see full materials list!" << G4endl;
277  return false;
278  }
279 
280  return true;
281 }
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
def SetMaterial
Definition: EmPlot.py:25
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetPhantomPosition ( G4ThreeVector  pos)

Definition at line 308 of file HadrontherapyDetectorConstruction.cc.

309 {
310  phantomPosition = pos;
311 }
static const G4double pos

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetPhantomSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 283 of file HadrontherapyDetectorConstruction.cc.

284 {
285  if (sizeX > 0.) phantomSizeX = sizeX;
286  if (sizeY > 0.) phantomSizeY = sizeY;
287  if (sizeZ > 0.) phantomSizeZ = sizeZ;
288 }

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::SetVoxelSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 300 of file HadrontherapyDetectorConstruction.cc.

301 {
302  if (sizeX > 0.) {sizeOfVoxelAlongX = sizeX;}
303  if (sizeY > 0.) {sizeOfVoxelAlongY = sizeY;}
304  if (sizeZ > 0.) {sizeOfVoxelAlongZ = sizeZ;}
305 }

Here is the caller graph for this function:

void HadrontherapyDetectorConstruction::UpdateGeometry ( )

Definition at line 320 of file HadrontherapyDetectorConstruction.cc.

321 {
322  /*
323  * Check parameters consistency
324  */
325  ParametersCheck();
326 
327  G4GeometryManager::GetInstance() -> OpenGeometry();
328  if (phantom)
329  {
330  phantom -> SetXHalfLength(phantomSizeX/2);
331  phantom -> SetYHalfLength(phantomSizeY/2);
332  phantom -> SetZHalfLength(phantomSizeZ/2);
333  phantomPhysicalVolume -> SetTranslation(phantomPosition);
334  }
335  else ConstructPhantom();
336 
337  // Get the center of the detector
339  if (detector)
340  {
341  detector -> SetXHalfLength(detectorSizeX/2);
342  detector -> SetYHalfLength(detectorSizeY/2);
343  detector -> SetZHalfLength(detectorSizeZ/2);
344  detectorPhysicalVolume -> SetTranslation(detectorPosition);
345  }
346  else ConstructDetector();
347 
348  // Round to nearest integer number of voxel
349 
350  numberOfVoxelsAlongX = G4lrint(detectorSizeX / sizeOfVoxelAlongX);
351  sizeOfVoxelAlongX = ( detectorSizeX / numberOfVoxelsAlongX );
352  numberOfVoxelsAlongY = G4lrint(detectorSizeY / sizeOfVoxelAlongY);
353  sizeOfVoxelAlongY = ( detectorSizeY / numberOfVoxelsAlongY );
354  numberOfVoxelsAlongZ = G4lrint(detectorSizeZ / sizeOfVoxelAlongZ);
355  sizeOfVoxelAlongZ = ( detectorSizeZ / numberOfVoxelsAlongZ );
357 
359 
361 
362  //Set parameters, either for the Construct() or for the UpdateROGeometry()
364  detectorSizeX/2,
365  detectorSizeY/2,
366  detectorSizeZ/2,
367  numberOfVoxelsAlongX,
368  numberOfVoxelsAlongY,
369  numberOfVoxelsAlongZ);
370 
371  //This method below has an effect only if the RO geometry is already built.
372  RO->UpdateROGeometry();
373 
374 
375 
376  volumeOfVoxel = sizeOfVoxelAlongX * sizeOfVoxelAlongY * sizeOfVoxelAlongZ;
377  massOfVoxel = detectorMaterial -> GetDensity() * volumeOfVoxel;
378  // This will clear the existing matrix (together with all data inside it)!
379  matrix = HadrontherapyMatrix::GetInstance(numberOfVoxelsAlongX,
380  numberOfVoxelsAlongY,
381  numberOfVoxelsAlongZ,
382  massOfVoxel);
383 
384 
385  // Comment out the line below if let calculation is not needed!
386  // Initialize LET with energy of primaries and clear data inside
387  if ( (let = HadrontherapyLet::GetInstance(this)) )
388  {
390  }
391 
392 
393  // Initialize analysis
394 #ifdef G4ANALYSIS_USE_ROOT
396  analysis -> flush(); // Finalize the root file
397  analysis -> book();
398 #endif
399  // Inform the kernel about the new geometry
401  G4RunManager::GetRunManager() -> PhysicsHasBeenModified();
402 
403  PrintParameters();
404  // CheckOverlaps();
405 }
static HadrontherapyAnalysisManager * GetInstance()
static HadrontherapyLet * GetInstance()
const G4VUserDetectorConstruction * GetUserDetectorConstruction() const
void Initialize()
Definition: errprop.cc:101
void Initialize(G4ThreeVector detectorPos, G4double detectorDimX, G4double detectorDimY, G4double detectorDimZ, G4int numberOfVoxelsX, G4int numberOfVoxelsY, G4int numberOfVoxelsZ)
static HadrontherapyMatrix * GetInstance()
G4VUserParallelWorld * GetParallelWorld(G4int i) const
static G4GeometryManager * GetInstance()
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

HadrontherapyDetectorSD* HadrontherapyDetectorConstruction::detectorSD

Definition at line 61 of file HadrontherapyDetectorConstruction.hh.

G4VPhysicalVolume* HadrontherapyDetectorConstruction::motherPhys

Definition at line 60 of file HadrontherapyDetectorConstruction.hh.


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