Geant4  10.02.p03
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 HadrontherapyDetectorConstructionGetInstance ()
 

Public Attributes

G4VPhysicalVolumemotherPhys
 
HadrontherapyDetectorSDdetectorSD
 

Private Member Functions

void ConstructPhantom ()
 
void ConstructDetector ()
 
void ParametersCheck ()
 
void CheckOverlaps ()
 

Private Attributes

HadrontherapyDetectorMessengerdetectorMessenger
 
G4VisAttributesskyBlue
 
G4VisAttributesred
 
HadrontherapyDetectorROGeometrydetectorROGeometry
 
HadrontherapyMatrixmatrix
 
HadrontherapyLetlet
 
G4Boxphantom
 
G4Boxdetector
 
G4LogicalVolumephantomLogicalVolume
 
G4LogicalVolumedetectorLogicalVolume
 
G4VPhysicalVolumephantomPhysicalVolume
 
G4VPhysicalVolumedetectorPhysicalVolume
 
G4double phantomSizeX
 
G4double phantomSizeY
 
G4double phantomSizeZ
 
G4double detectorSizeX
 
G4double detectorSizeY
 
G4double detectorSizeZ
 
G4ThreeVector phantomPosition
 
G4ThreeVector detectorPosition
 
G4ThreeVector detectorToPhantomPosition
 
G4double sizeOfVoxelAlongX
 
G4double sizeOfVoxelAlongY
 
G4double sizeOfVoxelAlongZ
 
G4int numberOfVoxelsAlongX
 
G4int numberOfVoxelsAlongY
 
G4int numberOfVoxelsAlongZ
 
G4double volumeOfVoxel
 
G4double massOfVoxel
 
G4MaterialphantomMaterial
 
G4MaterialdetectorMaterial
 
G4RegionaRegion
 

Static Private Attributes

static HadrontherapyDetectorConstructioninstance = 0
 

Detailed Description

Definition at line 48 of file HadrontherapyDetectorConstruction.hh.

Constructor & Destructor Documentation

◆ HadrontherapyDetectorConstruction()

HadrontherapyDetectorConstruction::HadrontherapyDetectorConstruction ( G4VPhysicalVolume physicalTreatmentRoom)

Definition at line 62 of file HadrontherapyDetectorConstruction.cc.

63 : motherPhys(physicalTreatmentRoom), // pointer to WORLD volume
65 phantom(0), detector(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
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 const double cm
Definition: G4SIunits.hh:118
static HadrontherapyAnalysisManager * GetInstance()
void SetPhantomSize(G4double sizeX, G4double sizeY, G4double sizeZ)
HadrontherapyDetectorMessenger * detectorMessenger
CLHEP::Hep3Vector G4ThreeVector
void SetDetectorSize(G4double sizeX, G4double sizeY, G4double sizeZ)
void SetDetectorToPhantomPosition(G4ThreeVector DetectorToPhantomPosition)
HadrontherapyDetectorROGeometry * detectorROGeometry
static const double um
Definition: G4SIunits.hh:112
Here is the call graph for this function:

◆ ~HadrontherapyDetectorConstruction()

HadrontherapyDetectorConstruction::~HadrontherapyDetectorConstruction ( )

Definition at line 102 of file HadrontherapyDetectorConstruction.cc.

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

Member Function Documentation

◆ CheckOverlaps()

void HadrontherapyDetectorConstruction::CheckOverlaps ( )
private

Definition at line 410 of file HadrontherapyDetectorConstruction.cc.

411 {
413  G4cout << thePVStore->size() << " physical volumes are defined" << G4endl;
414  G4bool overlapFlag = false;
415  G4int res=1000;
416  G4double tol=0.; //tolerance
417  for (size_t i=0;i<thePVStore->size();i++)
418  {
419  //overlapFlag = (*thePVStore)[i]->CheckOverlaps(res,tol,fCheckOverlapsVerbosity) | overlapFlag;
420  overlapFlag = (*thePVStore)[i]->CheckOverlaps(res,tol,true) | overlapFlag; }
421  if (overlapFlag)
422  G4cout << "Check: there are overlapping volumes" << G4endl;
423 }
int G4int
Definition: G4Types.hh:78
static G4PhysicalVolumeStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ConstructDetector()

void HadrontherapyDetectorConstruction::ConstructDetector ( )
private

Definition at line 170 of file HadrontherapyDetectorConstruction.cc.

171 {
172  // Definition of the solid volume of the Detector
173  detector = new G4Box("Detector",
174  detectorSizeX/2,
175  detectorSizeY/2,
176  detectorSizeZ/2);
177 
178  // Definition of the logic volume of the Phantom
181  "DetectorLog",
182  0,0,0);
183  // Definition of the physical volume of the Phantom
185  detectorPosition, // Setted by displacement
186  "DetectorPhys",
189  false,0);
190 
191  // Visualisation attributes of the detector
192  skyBlue = new G4VisAttributes( G4Colour(135/255. , 206/255. , 235/255. ));
193  skyBlue -> SetVisibility(true);
194  skyBlue -> SetForceSolid(true);
195  //skyBlue -> SetForceWireframe(true);
196  detectorLogicalVolume -> SetVisAttributes(skyBlue);
197 
198  // **************
199  // Cut per Region
200  // **************
201  //
202  // A smaller cut is fixed in the phantom to calculate the energy deposit with the
203  // required accuracy
204  if (!aRegion)
205  {
206  aRegion = new G4Region("DetectorLog");
207  detectorLogicalVolume -> SetRegion(aRegion);
208  aRegion -> AddRootLogicalVolume(detectorLogicalVolume);
209  }
210 }
Definition: G4Box.hh:64
Here is the caller graph for this function:

◆ ConstructPhantom()

void HadrontherapyDetectorConstruction::ConstructPhantom ( )
private

Definition at line 119 of file HadrontherapyDetectorConstruction.cc.

120 {
121  // Definition of the solid volume of the Phantom
122  phantom = new G4Box("Phantom",
123  phantomSizeX/2,
124  phantomSizeY/2,
125  phantomSizeZ/2);
126 
127  // Definition of the logical volume of the Phantom
130  "phantomLog", 0, 0, 0);
131 
132  // Definition of the physics volume of the Phantom
135  "phantomPhys",
137  motherPhys,
138  false,
139  0);
140 
141  // Visualisation attributes of the phantom
142  red = new G4VisAttributes(G4Colour(255/255., 0/255. ,0/255.));
143  red -> SetVisibility(true);
144  red -> SetForceSolid(true);
145  red -> SetForceWireframe(true);
146  phantomLogicalVolume -> SetVisAttributes(red);
147 }
Definition: G4Box.hh:64
Here is the caller graph for this function:

◆ GetDetectorLogicalVolume()

G4LogicalVolume* HadrontherapyDetectorConstruction::GetDetectorLogicalVolume ( )
inline

◆ GetDetectorToPhantomPosition()

G4ThreeVector HadrontherapyDetectorConstruction::GetDetectorToPhantomPosition ( )
inline

Definition at line 78 of file HadrontherapyDetectorConstruction.hh.

Here is the call graph for this function:

◆ GetDetectorToWorldPosition()

G4ThreeVector HadrontherapyDetectorConstruction::GetDetectorToWorldPosition ( )
inline

Definition at line 72 of file HadrontherapyDetectorConstruction.hh.

Here is the caller graph for this function:

◆ GetInstance()

HadrontherapyDetectorConstruction * HadrontherapyDetectorConstruction::GetInstance ( )
static

Definition at line 110 of file HadrontherapyDetectorConstruction.cc.

111 {
112  return instance;
113 }
static HadrontherapyDetectorConstruction * instance

◆ InitializeDetectorROGeometry()

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,
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:

◆ IsInside()

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
G4GLOB_DLL std::ostream G4cout
double getX() const
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:

◆ ParametersCheck()

void HadrontherapyDetectorConstruction::ParametersCheck ( )
private

Definition at line 228 of file HadrontherapyDetectorConstruction.cc.

229 {
230  // Check phantom/detector sizes & relative position
231  if (!IsInside(detectorSizeX,
234  phantomSizeX,
235  phantomSizeY,
236  phantomSizeZ,
238  ))
239  G4Exception("HadrontherapyDetectorConstruction::ParametersCheck()", "Hadrontherapy0001", FatalException, "Error: Detector is not fully inside Phantom!");
240 
241  // Check Detector sizes respect to the voxel ones
242 
244  G4Exception("HadrontherapyDetectorConstruction::ParametersCheck()", "Hadrontherapy0002", FatalException, "Error: Detector X size must be bigger or equal than that of Voxel X!");
245  }
247  G4Exception(" HadrontherapyDetectorConstruction::ParametersCheck()", "Hadrontherapy0003", FatalException, "Error: Detector Y size must be bigger or equal than that of Voxel Y!");
248  }
250  G4Exception(" HadrontherapyDetectorConstruction::ParametersCheck()", "Hadrontherapy0004", FatalException, "Error: Detector Z size must be bigger or equal than that of Voxel Z!");
251  }
252 }
bool IsInside(G4double detectorX, G4double detectorY, G4double detectorZ, G4double phantomX, G4double phantomY, G4double phantomZ, G4ThreeVector pos)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PrintParameters()

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 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
double getY() const
G4GLOB_DLL std::ostream G4cout
double getX() const
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:

◆ SetDetectorPosition()

void HadrontherapyDetectorConstruction::SetDetectorPosition ( )
inline

Definition at line 88 of file HadrontherapyDetectorConstruction.hh.

89  {
90  // Adjust detector position
94 
95 
96  }
void setY(double)
void setZ(double)
void setX(double)
double getY() const
double getX() const
double getZ() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetDetectorSize()

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

Definition at line 291 of file HadrontherapyDetectorConstruction.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetDetectorToPhantomPosition()

void HadrontherapyDetectorConstruction::SetDetectorToPhantomPosition ( G4ThreeVector  DetectorToPhantomPosition)

Definition at line 314 of file HadrontherapyDetectorConstruction.cc.

Here is the caller graph for this function:

◆ SetPhantomMaterial()

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;
263  {
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
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
#define G4endl
Definition: G4ios.hh:61
def SetMaterial(material_name)
Definition: EmPlot.py:25
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPhantomPosition()

void HadrontherapyDetectorConstruction::SetPhantomPosition ( G4ThreeVector  pos)

Definition at line 308 of file HadrontherapyDetectorConstruction.cc.

309 {
311 }
static const G4double pos
Here is the caller graph for this function:

◆ SetPhantomSize()

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:

◆ SetVoxelSize()

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:

◆ UpdateGeometry()

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 
357 
359 
361 
362  //Set parameters, either for the Construct() or for the UpdateROGeometry()
364  detectorSizeX/2,
365  detectorSizeY/2,
366  detectorSizeZ/2,
370 
371  //This method below has an effect only if the RO geometry is already built.
372  RO->UpdateROGeometry();
373 
374 
375 
377  massOfVoxel = detectorMaterial -> GetDensity() * volumeOfVoxel;
378  // This will clear the existing matrix (together with all data inside it)!
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()
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()
static G4GeometryManager * GetInstance()
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
int G4lrint(double ad)
Definition: templates.hh:163
G4VUserParallelWorld * GetParallelWorld(G4int i) const
const G4VUserDetectorConstruction * GetUserDetectorConstruction() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ aRegion

G4Region* HadrontherapyDetectorConstruction::aRegion
private

Definition at line 198 of file HadrontherapyDetectorConstruction.hh.

◆ detector

G4Box * HadrontherapyDetectorConstruction::detector
private

Definition at line 173 of file HadrontherapyDetectorConstruction.hh.

◆ detectorLogicalVolume

G4LogicalVolume * HadrontherapyDetectorConstruction::detectorLogicalVolume
private

Definition at line 174 of file HadrontherapyDetectorConstruction.hh.

◆ detectorMaterial

G4Material * HadrontherapyDetectorConstruction::detectorMaterial
private

Definition at line 197 of file HadrontherapyDetectorConstruction.hh.

◆ detectorMessenger

HadrontherapyDetectorMessenger* HadrontherapyDetectorConstruction::detectorMessenger
private

Definition at line 164 of file HadrontherapyDetectorConstruction.hh.

◆ detectorPhysicalVolume

G4VPhysicalVolume * HadrontherapyDetectorConstruction::detectorPhysicalVolume
private

Definition at line 175 of file HadrontherapyDetectorConstruction.hh.

◆ detectorPosition

G4ThreeVector HadrontherapyDetectorConstruction::detectorPosition
private

Definition at line 185 of file HadrontherapyDetectorConstruction.hh.

◆ detectorROGeometry

HadrontherapyDetectorROGeometry* HadrontherapyDetectorConstruction::detectorROGeometry
private

Definition at line 169 of file HadrontherapyDetectorConstruction.hh.

◆ detectorSD

HadrontherapyDetectorSD* HadrontherapyDetectorConstruction::detectorSD

Definition at line 61 of file HadrontherapyDetectorConstruction.hh.

◆ detectorSizeX

G4double HadrontherapyDetectorConstruction::detectorSizeX
private

Definition at line 181 of file HadrontherapyDetectorConstruction.hh.

◆ detectorSizeY

G4double HadrontherapyDetectorConstruction::detectorSizeY
private

Definition at line 182 of file HadrontherapyDetectorConstruction.hh.

◆ detectorSizeZ

G4double HadrontherapyDetectorConstruction::detectorSizeZ
private

Definition at line 183 of file HadrontherapyDetectorConstruction.hh.

◆ detectorToPhantomPosition

G4ThreeVector HadrontherapyDetectorConstruction::detectorToPhantomPosition
private

Definition at line 185 of file HadrontherapyDetectorConstruction.hh.

◆ instance

HadrontherapyDetectorConstruction * HadrontherapyDetectorConstruction::instance = 0
staticprivate

Definition at line 163 of file HadrontherapyDetectorConstruction.hh.

◆ let

HadrontherapyLet* HadrontherapyDetectorConstruction::let
private

Definition at line 171 of file HadrontherapyDetectorConstruction.hh.

◆ massOfVoxel

G4double HadrontherapyDetectorConstruction::massOfVoxel
private

Definition at line 195 of file HadrontherapyDetectorConstruction.hh.

◆ matrix

HadrontherapyMatrix* HadrontherapyDetectorConstruction::matrix
private

Definition at line 170 of file HadrontherapyDetectorConstruction.hh.

◆ motherPhys

G4VPhysicalVolume* HadrontherapyDetectorConstruction::motherPhys

Definition at line 60 of file HadrontherapyDetectorConstruction.hh.

◆ numberOfVoxelsAlongX

G4int HadrontherapyDetectorConstruction::numberOfVoxelsAlongX
private

Definition at line 191 of file HadrontherapyDetectorConstruction.hh.

◆ numberOfVoxelsAlongY

G4int HadrontherapyDetectorConstruction::numberOfVoxelsAlongY
private

Definition at line 192 of file HadrontherapyDetectorConstruction.hh.

◆ numberOfVoxelsAlongZ

G4int HadrontherapyDetectorConstruction::numberOfVoxelsAlongZ
private

Definition at line 193 of file HadrontherapyDetectorConstruction.hh.

◆ phantom

G4Box* HadrontherapyDetectorConstruction::phantom
private

Definition at line 173 of file HadrontherapyDetectorConstruction.hh.

◆ phantomLogicalVolume

G4LogicalVolume* HadrontherapyDetectorConstruction::phantomLogicalVolume
private

Definition at line 174 of file HadrontherapyDetectorConstruction.hh.

◆ phantomMaterial

G4Material* HadrontherapyDetectorConstruction::phantomMaterial
private

Definition at line 197 of file HadrontherapyDetectorConstruction.hh.

◆ phantomPhysicalVolume

G4VPhysicalVolume* HadrontherapyDetectorConstruction::phantomPhysicalVolume
private

Definition at line 175 of file HadrontherapyDetectorConstruction.hh.

◆ phantomPosition

G4ThreeVector HadrontherapyDetectorConstruction::phantomPosition
private

Definition at line 185 of file HadrontherapyDetectorConstruction.hh.

◆ phantomSizeX

G4double HadrontherapyDetectorConstruction::phantomSizeX
private

Definition at line 177 of file HadrontherapyDetectorConstruction.hh.

◆ phantomSizeY

G4double HadrontherapyDetectorConstruction::phantomSizeY
private

Definition at line 178 of file HadrontherapyDetectorConstruction.hh.

◆ phantomSizeZ

G4double HadrontherapyDetectorConstruction::phantomSizeZ
private

Definition at line 179 of file HadrontherapyDetectorConstruction.hh.

◆ red

G4VisAttributes* HadrontherapyDetectorConstruction::red
private

Definition at line 167 of file HadrontherapyDetectorConstruction.hh.

◆ sizeOfVoxelAlongX

G4double HadrontherapyDetectorConstruction::sizeOfVoxelAlongX
private

Definition at line 187 of file HadrontherapyDetectorConstruction.hh.

◆ sizeOfVoxelAlongY

G4double HadrontherapyDetectorConstruction::sizeOfVoxelAlongY
private

Definition at line 188 of file HadrontherapyDetectorConstruction.hh.

◆ sizeOfVoxelAlongZ

G4double HadrontherapyDetectorConstruction::sizeOfVoxelAlongZ
private

Definition at line 189 of file HadrontherapyDetectorConstruction.hh.

◆ skyBlue

G4VisAttributes* HadrontherapyDetectorConstruction::skyBlue
private

Definition at line 166 of file HadrontherapyDetectorConstruction.hh.

◆ volumeOfVoxel

G4double HadrontherapyDetectorConstruction::volumeOfVoxel
private

Definition at line 195 of file HadrontherapyDetectorConstruction.hh.


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