Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions | Data Fields
HadrontherapyDetectorConstruction Class Reference

#include <HadrontherapyDetectorConstruction.hh>

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 ()
 

Data Fields

G4VPhysicalVolumemotherPhys
 
HadrontherapyDetectorSDdetectorSD
 

Detailed Description

Definition at line 52 of file HadrontherapyDetectorConstruction.hh.

Constructor & Destructor Documentation

HadrontherapyDetectorConstruction::HadrontherapyDetectorConstruction ( G4VPhysicalVolume physicalTreatmentRoom)

Definition at line 70 of file HadrontherapyDetectorConstruction.cc.

References python.hepunit::cm, HadrontherapyAnalysisManager::GetInstance(), SetDetectorPosition(), SetDetectorSize(), SetDetectorToPhantomPosition(), SetPhantomMaterial(), SetPhantomPosition(), SetPhantomSize(), and UpdateGeometry().

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

Definition at line 111 of file HadrontherapyDetectorConstruction.cc.

112 {
113  delete detectorROGeometry;
114  delete matrix;
115  delete detectorMessenger;
116 }

Member Function Documentation

G4LogicalVolume* HadrontherapyDetectorConstruction::GetDetectorLogicalVolume ( )
inline

Definition at line 163 of file HadrontherapyDetectorConstruction.hh.

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

Definition at line 81 of file HadrontherapyDetectorConstruction.hh.

References CLHEP::Hep3Vector::getX(), CLHEP::Hep3Vector::getY(), and CLHEP::Hep3Vector::getZ().

82 {
83  return G4ThreeVector(phantomSizeX/2 - detectorSizeX/2 + detectorPosition.getX(),
84  phantomSizeY/2 - detectorSizeY/2 + detectorPosition.getY(),
85  phantomSizeZ/2 - detectorSizeZ/2 + detectorPosition.getZ()
86  );
87 }
CLHEP::Hep3Vector G4ThreeVector
double getY() const
double getX() const
double getZ() const
G4ThreeVector HadrontherapyDetectorConstruction::GetDetectorToWorldPosition ( )
inline

Definition at line 75 of file HadrontherapyDetectorConstruction.hh.

Referenced by PassiveCarbonBeamLine::Construct(), PassiveProtonBeamLine::Construct(), and UpdateGeometry().

76  {
77  return phantomPosition + detectorPosition;
78  }
HadrontherapyDetectorConstruction * HadrontherapyDetectorConstruction::GetInstance ( )
static

Definition at line 119 of file HadrontherapyDetectorConstruction.cc.

120 {
121 
122  return instance;
123 }
void HadrontherapyDetectorConstruction::InitializeDetectorROGeometry ( HadrontherapyDetectorROGeometry RO,
G4ThreeVector  detectorToWorldPosition 
)

Definition at line 230 of file HadrontherapyDetectorConstruction.cc.

References HadrontherapyDetectorROGeometry::Initialize().

Referenced by PassiveCarbonBeamLine::Construct(), and PassiveProtonBeamLine::Construct().

234 {
235  RO->Initialize(detectorToWorldPosition,
236  detectorSizeX/2,
237  detectorSizeY/2,
238  detectorSizeZ/2,
239  numberOfVoxelsAlongX,
240  numberOfVoxelsAlongY,
241  numberOfVoxelsAlongZ);
242 }
void Initialize(G4ThreeVector detectorPos, G4double detectorDimX, G4double detectorDimY, G4double detectorDimZ, G4int numberOfVoxelsX, G4int numberOfVoxelsY, G4int numberOfVoxelsZ)
bool HadrontherapyDetectorConstruction::IsInside ( G4double  detectorX,
G4double  detectorY,
G4double  detectorZ,
G4double  phantomX,
G4double  phantomY,
G4double  phantomZ,
G4ThreeVector  pos 
)
inline

Definition at line 102 of file HadrontherapyDetectorConstruction.hh.

References G4cout, G4endl, CLHEP::Hep3Vector::getX(), CLHEP::Hep3Vector::getY(), and CLHEP::Hep3Vector::getZ().

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

Definition at line 433 of file HadrontherapyDetectorConstruction.cc.

References G4BestUnit, G4cout, G4endl, CLHEP::Hep3Vector::getX(), CLHEP::Hep3Vector::getY(), and CLHEP::Hep3Vector::getZ().

Referenced by UpdateGeometry().

434 {
435 
436  G4cout << "The (X,Y,Z) dimensions of the phantom are : (" <<
437  G4BestUnit( phantom -> GetXHalfLength()*2., "Length") << ',' <<
438  G4BestUnit( phantom -> GetYHalfLength()*2., "Length") << ',' <<
439  G4BestUnit( phantom -> GetZHalfLength()*2., "Length") << ')' << G4endl;
440 
441  G4cout << "The (X,Y,Z) dimensions of the detector are : (" <<
442  G4BestUnit( detector -> GetXHalfLength()*2., "Length") << ',' <<
443  G4BestUnit( detector -> GetYHalfLength()*2., "Length") << ',' <<
444  G4BestUnit( detector -> GetZHalfLength()*2., "Length") << ')' << G4endl;
445 
446  G4cout << "Displacement between Phantom and World is: ";
447  G4cout << "DX= "<< G4BestUnit(phantomPosition.getX(),"Length") <<
448  "DY= "<< G4BestUnit(phantomPosition.getY(),"Length") <<
449  "DZ= "<< G4BestUnit(phantomPosition.getZ(),"Length") << G4endl;
450 
451  G4cout << "The (X,Y,Z) sizes of the Voxels are: (" <<
452  G4BestUnit(sizeOfVoxelAlongX, "Length") << ',' <<
453  G4BestUnit(sizeOfVoxelAlongY, "Length") << ',' <<
454  G4BestUnit(sizeOfVoxelAlongZ, "Length") << ')' << G4endl;
455 
456  G4cout << "The number of Voxels along (X,Y,Z) is: (" <<
457  numberOfVoxelsAlongX << ',' <<
458  numberOfVoxelsAlongY <<',' <<
459  numberOfVoxelsAlongZ << ')' << G4endl;
460 
461 }
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
void HadrontherapyDetectorConstruction::SetDetectorPosition ( )
inline

Definition at line 91 of file HadrontherapyDetectorConstruction.hh.

References CLHEP::Hep3Vector::getX(), CLHEP::Hep3Vector::getY(), CLHEP::Hep3Vector::getZ(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), and CLHEP::Hep3Vector::setZ().

Referenced by HadrontherapyDetectorConstruction(), and UpdateGeometry().

92  {
93  // Adjust detector position
94  detectorPosition.setX(detectorToPhantomPosition.getX() - phantomSizeX/2 + detectorSizeX/2);
95  detectorPosition.setY(detectorToPhantomPosition.getY() - phantomSizeY/2 + detectorSizeY/2);
96  detectorPosition.setZ(detectorToPhantomPosition.getZ() - phantomSizeZ/2 + detectorSizeZ/2);
97 
98 
99  }
double getY() const
void setY(double)
void setZ(double)
void setX(double)
double getX() const
double getZ() const
void HadrontherapyDetectorConstruction::SetDetectorSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 316 of file HadrontherapyDetectorConstruction.cc.

References SetVoxelSize().

Referenced by HadrontherapyDetectorConstruction().

317 {
318  if (sizeX > 0.) {detectorSizeX = sizeX;}
319  if (sizeY > 0.) {detectorSizeY = sizeY;}
320  if (sizeZ > 0.) {detectorSizeZ = sizeZ;}
321  SetVoxelSize(sizeOfVoxelAlongX, sizeOfVoxelAlongY, sizeOfVoxelAlongZ);
322 }
void SetVoxelSize(G4double sizeX, G4double sizeY, G4double sizeZ)
void HadrontherapyDetectorConstruction::SetDetectorToPhantomPosition ( G4ThreeVector  DetectorToPhantomPosition)

Definition at line 337 of file HadrontherapyDetectorConstruction.cc.

Referenced by HadrontherapyDetectorConstruction().

338 {
339  detectorToPhantomPosition = displ;
340 }
G4bool HadrontherapyDetectorConstruction::SetPhantomMaterial ( G4String  material)

Definition at line 280 of file HadrontherapyDetectorConstruction.cc.

References G4cout, G4endl, G4RunManager::GetRunManager(), G4NistManager::Instance(), and EmPlot::SetMaterial().

Referenced by HadrontherapyDetectorConstruction().

281 {
282 
283  if (G4Material* pMat = G4NistManager::Instance()->FindOrBuildMaterial(material, false) )
284  {
285  phantomMaterial = pMat;
286  detectorMaterial = pMat;
287  if (detectorLogicalVolume && phantomLogicalVolume)
288  {
289  detectorLogicalVolume -> SetMaterial(pMat);
290  phantomLogicalVolume -> SetMaterial(pMat);
291 
292  G4RunManager::GetRunManager() -> PhysicsHasBeenModified();
293  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
294  G4cout << "The material of Phantom/Detector has been changed to " << material << G4endl;
295  }
296  }
297  else
298  {
299  G4cout << "WARNING: material \"" << material << "\" doesn't exist in NIST elements/materials"
300  " table [located in $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl;
301  G4cout << "Use command \"/parameter/nist\" to see full materials list!" << G4endl;
302  return false;
303  }
304 
305  return true;
306 }
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
def SetMaterial
Definition: EmPlot.py:25
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
#define G4endl
Definition: G4ios.hh:61
void HadrontherapyDetectorConstruction::SetPhantomPosition ( G4ThreeVector  pos)

Definition at line 331 of file HadrontherapyDetectorConstruction.cc.

Referenced by HadrontherapyDetectorConstruction().

332 {
333  phantomPosition = pos;
334 }
void HadrontherapyDetectorConstruction::SetPhantomSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 308 of file HadrontherapyDetectorConstruction.cc.

Referenced by HadrontherapyDetectorConstruction().

309 {
310  if (sizeX > 0.) phantomSizeX = sizeX;
311  if (sizeY > 0.) phantomSizeY = sizeY;
312  if (sizeZ > 0.) phantomSizeZ = sizeZ;
313 }
void HadrontherapyDetectorConstruction::SetVoxelSize ( G4double  sizeX,
G4double  sizeY,
G4double  sizeZ 
)

Definition at line 325 of file HadrontherapyDetectorConstruction.cc.

Referenced by SetDetectorSize().

326 {
327  if (sizeX > 0.) {sizeOfVoxelAlongX = sizeX;}
328  if (sizeY > 0.) {sizeOfVoxelAlongY = sizeY;}
329  if (sizeZ > 0.) {sizeOfVoxelAlongZ = sizeZ;}
330 }
void HadrontherapyDetectorConstruction::UpdateGeometry ( )

Definition at line 342 of file HadrontherapyDetectorConstruction.cc.

References G4lrint(), GetDetectorToWorldPosition(), HadrontherapyLet::GetInstance(), HadrontherapyAnalysisManager::GetInstance(), HadrontherapyMatrix::GetInstance(), G4GeometryManager::GetInstance(), G4VUserDetectorConstruction::GetParallelWorld(), G4RunManager::GetRunManager(), G4RunManager::GetUserDetectorConstruction(), HadrontherapyDetectorROGeometry::Initialize(), Initialize(), PrintParameters(), SetDetectorPosition(), and HadrontherapyDetectorROGeometry::UpdateROGeometry().

Referenced by HadrontherapyDetectorConstruction().

343 {
344  /*
345  * Check parameters consistency
346  */
347  ParametersCheck();
348 
349  G4GeometryManager::GetInstance() -> OpenGeometry();
350  if (phantom)
351  {
352  phantom -> SetXHalfLength(phantomSizeX/2);
353  phantom -> SetYHalfLength(phantomSizeY/2);
354  phantom -> SetZHalfLength(phantomSizeZ/2);
355  phantomPhysicalVolume -> SetTranslation(phantomPosition);
356  }
357  else ConstructPhantom();
358 
359  // Get the center of the detector
361  if (detector)
362  {
363  detector -> SetXHalfLength(detectorSizeX/2);
364  detector -> SetYHalfLength(detectorSizeY/2);
365  detector -> SetZHalfLength(detectorSizeZ/2);
366  detectorPhysicalVolume -> SetTranslation(detectorPosition);
367  }
368  else ConstructDetector();
369 
370  // Round to nearest integer number of voxel
371 
372  numberOfVoxelsAlongX = G4lrint(detectorSizeX / sizeOfVoxelAlongX);
373 
374  sizeOfVoxelAlongX = ( detectorSizeX / numberOfVoxelsAlongX );
375 
376  numberOfVoxelsAlongY = G4lrint(detectorSizeY / sizeOfVoxelAlongY);
377 
378  sizeOfVoxelAlongY = ( detectorSizeY / numberOfVoxelsAlongY );
379 
380  numberOfVoxelsAlongZ = G4lrint(detectorSizeZ / sizeOfVoxelAlongZ);
381 
382  sizeOfVoxelAlongZ = ( detectorSizeZ / numberOfVoxelsAlongZ );
383 
385 
387 
388 
390 
391  //Set parameters, either for the Construct() or for the UpdateROGeometry()
393  detectorSizeX/2,
394  detectorSizeY/2,
395  detectorSizeZ/2,
396  numberOfVoxelsAlongX,
397  numberOfVoxelsAlongY,
398  numberOfVoxelsAlongZ);
399 
400  //This method below has an effect only if the RO geometry is already built.
401  RO->UpdateROGeometry();
402 
403  volumeOfVoxel = sizeOfVoxelAlongX * sizeOfVoxelAlongY * sizeOfVoxelAlongZ;
404  massOfVoxel = detectorMaterial -> GetDensity() * volumeOfVoxel;
405  // This will clear the existing matrix (together with all data inside it)!
406  matrix = HadrontherapyMatrix::GetInstance(numberOfVoxelsAlongX,
407  numberOfVoxelsAlongY,
408  numberOfVoxelsAlongZ,
409  massOfVoxel);
410 
411 
412 // Comment out the line below if let calculation is not needed!
413  // Initialize LET with energy of primaries and clear data inside
414  if ( (let = HadrontherapyLet::GetInstance(this)) )
415  {
417  }
418 
419 
420  // Initialize analysis
421 #ifdef G4ANALYSIS_USE_ROOT
423  analysis -> flush(); // Finalize the root file
424  analysis -> book();
425 #endif
426  // Inform the kernel about the new geometry
427  G4RunManager::GetRunManager() -> GeometryHasBeenModified();
428  G4RunManager::GetRunManager() -> PhysicsHasBeenModified();
429 
430  PrintParameters();
431 }
static HadrontherapyAnalysisManager * GetInstance()
static HadrontherapyLet * GetInstance()
const G4VUserDetectorConstruction * GetUserDetectorConstruction() const
void Initialize()
Definition: errprop.cc:96
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:74
int G4lrint(double ad)
Definition: templates.hh:163

Field Documentation

HadrontherapyDetectorSD* HadrontherapyDetectorConstruction::detectorSD

Definition at line 65 of file HadrontherapyDetectorConstruction.hh.

G4VPhysicalVolume* HadrontherapyDetectorConstruction::motherPhys

Definition at line 64 of file HadrontherapyDetectorConstruction.hh.


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