Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
DicomDetectorConstruction Class Referenceabstract

#include <DicomDetectorConstruction.hh>

Inheritance diagram for DicomDetectorConstruction:
G4VUserDetectorConstruction DicomNestedParamDetectorConstruction DicomPartialDetectorConstruction DicomRegularDetectorConstruction

Public Member Functions

 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
 

Protected Member Functions

void InitialisationOfMaterials ()
 
void ReadPhantomData ()
 
void ReadPhantomDataFile (const G4String &fname)
 
void MergeZSliceHeaders ()
 
G4MaterialBuildMaterialWithChangingDensity (const G4Material *origMate, float density, G4String newMateName)
 
void ConstructPhantomContainer ()
 
virtual void ConstructPhantom ()=0
 
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

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
 
DicomPhantomZSliceMergedmergedSlices
 
std::set< G4LogicalVolume * > scorers
 
G4bool fConstructed
 

Detailed Description

Dicom detector construction

 - Start the building of the geometry
 - Initialisation of materials
 - Creation of the world
 - Reading of the DICOM data

History: 30.11.07 First version

Author
P. Arce

Definition at line 58 of file DicomDetectorConstruction.hh.

Constructor & Destructor Documentation

DicomDetectorConstruction::DicomDetectorConstruction ( )

Definition at line 53 of file DicomDetectorConstruction.cc.

55  fAir(0),
56 
57  fWorld_solid(0),
58  fWorld_logic(0),
59  fWorld_phys(0),
60 
63  fContainer_phys(0),
64 
65  fNoFiles(0),
66  fMateIDs(0),
67 
69 
70  fNVoxelX(0),
71  fNVoxelY(0),
72  fNVoxelZ(0),
73  fVoxelHalfDimX(0),
74  fVoxelHalfDimY(0),
75  fVoxelHalfDimZ(0),
76 
77  fConstructed(false)
78 {
79 
80 }
DicomPhantomZSliceHeader * fZSliceHeaderMerged
DicomDetectorConstruction::~DicomDetectorConstruction ( )

Definition at line 83 of file DicomDetectorConstruction.cc.

84 {
85 }

Member Function Documentation

G4Material * DicomDetectorConstruction::BuildMaterialWithChangingDensity ( const G4Material origMate,
float  density,
G4String  newMateName 
)
protected

Definition at line 465 of file DicomDetectorConstruction.cc.

References G4Material::AddElement(), python.hepunit::cm3, g(), G4Material::GetElement(), G4Material::GetFractionVector(), G4Material::GetNumberOfElements(), kStateUndefined, and python.hepunit::STP_Temperature.

Referenced by ReadPhantomDataFile().

467 {
468  //----- Copy original material, but with new density
469  G4int nelem = origMate->GetNumberOfElements();
470  G4Material* mate = new G4Material( newMateName, density*g/cm3, nelem,
472 
473  for( G4int ii = 0; ii < nelem; ii++ ){
474  G4double frac = origMate->GetFractionVector()[ii];
475  G4Element* elem = const_cast<G4Element*>(origMate->GetElement(ii));
476  mate->AddElement( elem, frac );
477  }
478 
479  return mate;
480 }
float STP_Temperature
Definition: hepunit.py:302
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:345
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
const G4double * GetFractionVector() const
Definition: G4Material.hh:192
G4VPhysicalVolume * DicomDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Reimplemented in DicomPartialDetectorConstruction.

Definition at line 88 of file DicomDetectorConstruction.cc.

References ConstructPhantom(), ConstructPhantomContainer(), fAir, fConstructed, fWorld_logic, fWorld_phys, fWorld_solid, InitialisationOfMaterials(), python.hepunit::m, and ReadPhantomData().

89 {
90  if(!fConstructed || fWorld_phys == 0) {
91  fConstructed = true;
93 
94  //----- Build world
95  G4double worldXDimension = 1.*m;
96  G4double worldYDimension = 1.*m;
97  G4double worldZDimension = 1.*m;
98 
99  fWorld_solid = new G4Box( "WorldSolid",
100  worldXDimension,
101  worldYDimension,
102  worldZDimension );
103 
105  fAir,
106  "WorldLogical",
107  0, 0, 0 );
108 
109  fWorld_phys = new G4PVPlacement( 0,
110  G4ThreeVector(0,0,0),
111  "World",
112  fWorld_logic,
113  0,
114  false,
115  0 );
116 
117  ReadPhantomData();
118 
121  }
122  return fWorld_phys;
123 }
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:63
double G4double
Definition: G4Types.hh:76
virtual void ConstructPhantom()=0
virtual void DicomDetectorConstruction::ConstructPhantom ( )
protectedpure virtual

Referenced by Construct().

void DicomDetectorConstruction::ConstructPhantomContainer ( )
protected

Definition at line 483 of file DicomDetectorConstruction.cc.

References fContainer_logic, fContainer_phys, fContainer_solid, fMaterials, fNVoxelX, fNVoxelY, fNVoxelZ, fVoxelHalfDimX, fVoxelHalfDimY, fVoxelHalfDimZ, fWorld_logic, fZSliceHeaderMerged, G4cout, G4endl, DicomPhantomZSliceHeader::GetMaxX(), DicomPhantomZSliceHeader::GetMaxY(), DicomPhantomZSliceHeader::GetMaxZ(), DicomPhantomZSliceHeader::GetMinX(), DicomPhantomZSliceHeader::GetMinY(), DicomPhantomZSliceHeader::GetMinZ(), DicomPhantomZSliceHeader::GetNoVoxelX(), DicomPhantomZSliceHeader::GetNoVoxelY(), DicomPhantomZSliceHeader::GetNoVoxelZ(), DicomPhantomZSliceHeader::GetVoxelHalfX(), DicomPhantomZSliceHeader::GetVoxelHalfY(), and DicomPhantomZSliceHeader::GetVoxelHalfZ().

Referenced by Construct().

484 {
485  //---- Extract number of voxels and voxel dimensions
489 
493 #ifdef G4VERBOSE
494  G4cout << " fNVoxelX " << fNVoxelX << " fVoxelHalfDimX " << fVoxelHalfDimX <<G4endl;
495  G4cout << " fNVoxelY " << fNVoxelY << " fVoxelHalfDimY " << fVoxelHalfDimY <<G4endl;
496  G4cout << " fNVoxelZ " << fNVoxelZ << " fVoxelHalfDimZ " << fVoxelHalfDimZ <<G4endl;
497  G4cout << " totalPixels " << fNVoxelX*fNVoxelY*fNVoxelZ << G4endl;
498 #endif
499 
500  //----- Define the volume that contains all the voxels
501  fContainer_solid = new G4Box("phantomContainer",fNVoxelX*fVoxelHalfDimX,
506  //the material is not important, it will be fully filled by the voxels
507  fMaterials[0],
508  "phantomContainer",
509  0, 0, 0 );
510  //--- Place it on the world
514  G4ThreeVector posCentreVoxels(fOffsetX,fOffsetY,fOffsetZ);
515 #ifdef G4VERBOSE
516  G4cout << " placing voxel container volume at " << posCentreVoxels << G4endl;
517 #endif
519  new G4PVPlacement(0, // rotation
520  posCentreVoxels,
521  fContainer_logic, // The logic volume
522  "phantomContainer", // Name
523  fWorld_logic, // Mother
524  false, // No op. bool.
525  1); // Copy number
526 
527  //fContainer_logic->SetVisAttributes(new G4VisAttributes(G4Colour(1.,0.,0.)));
528 }
Definition: G4Box.hh:63
DicomPhantomZSliceHeader * fZSliceHeaderMerged
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
std::vector< G4Material * > fMaterials
void DicomDetectorConstruction::ConstructSDandField ( )
protectedvirtual

Reimplemented from G4VUserDetectorConstruction.

Definition at line 568 of file DicomDetectorConstruction.cc.

References fNVoxelX, fNVoxelY, fNVoxelZ, G4cout, G4endl, G4MultiFunctionalDetector::RegisterPrimitive(), scorers, and G4VUserDetectorConstruction::SetSensitiveDetector().

569 {
570 
571  G4cout << "\n\n\n\n\t CONSTRUCT SD AND FIELD \n\n\n" << G4endl;
572 
573  //G4SDManager* SDman = G4SDManager::GetSDMpointer();
574 
575  //SDman->SetVerboseLevel(1);
576 
577  //
578  // Sensitive Detector Name
579  G4String concreteSDname = "phantomSD";
580  std::vector<G4String> scorer_names;
581  scorer_names.push_back(concreteSDname);
582  //------------------------
583  // MultiFunctionalDetector
584  //------------------------
585  //
586  // Define MultiFunctionalDetector with name.
587  // declare MFDet as a MultiFunctionalDetector scorer
588  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
589  //SDman->AddNewDetector( MFDet ); // Register SD to SDManager
590  //G4VPrimitiveScorer* dosedep = new G4PSDoseDeposit("DoseDeposit");
591  G4VPrimitiveScorer* dosedep = new G4PSDoseDeposit3D("DoseDeposit", fNVoxelX, fNVoxelY,
592  fNVoxelZ);
593  MFDet->RegisterPrimitive(dosedep);
594 
595  for(std::set<G4LogicalVolume*>::iterator ite = scorers.begin(); ite != scorers.end(); ++ite) {
596  SetSensitiveDetector(*ite, MFDet);
597  }
598 
599  /*if(DicomRunAction::Instance()->GetDicomRun()) {
600  DicomRunAction::Instance()->GetDicomRun()->ConstructMFD(scorer_names);
601  }*/
602 
603 
604 }
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
std::set< G4LogicalVolume * > scorers
G4GLOB_DLL std::ostream G4cout
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
#define G4endl
Definition: G4ios.hh:61
void DicomDetectorConstruction::InitialisationOfMaterials ( )
protected

Definition at line 127 of file DicomDetectorConstruction.cc.

References test::a, G4Material::AddElement(), python.hepunit::cm3, density, elC, elH, elN, elO, fAir, fOriginalMaterials, g(), python.hepunit::mg, python.hepunit::mole, symbol, and z.

Referenced by DicomPartialDetectorConstruction::Construct(), and Construct().

128 {
129  // Creating elements :
130  G4double z, a, density;
132 
133  G4Element* elC = new G4Element( name = "Carbon",
134  symbol = "C",
135  z = 6.0, a = 12.011 * g/mole );
136  G4Element* elH = new G4Element( name = "Hydrogen",
137  symbol = "H",
138  z = 1.0, a = 1.008 * g/mole );
139  G4Element* elN = new G4Element( name = "Nitrogen",
140  symbol = "N",
141  z = 7.0, a = 14.007 * g/mole );
142  G4Element* elO = new G4Element( name = "Oxygen",
143  symbol = "O",
144  z = 8.0, a = 16.00 * g/mole );
145  G4Element* elNa = new G4Element( name = "Sodium",
146  symbol = "Na",
147  z= 11.0, a = 22.98977* g/mole );
148  G4Element* elS = new G4Element( name = "Sulfur",
149  symbol = "S",
150  z = 16.0,a = 32.065* g/mole );
151  G4Element* elCl = new G4Element( name = "Chlorine",
152  symbol = "P",
153  z = 17.0, a = 35.453* g/mole );
154  G4Element* elK = new G4Element( name = "Potassium",
155  symbol = "P",
156  z = 19.0, a = 30.0983* g/mole );
157  G4Element* elP = new G4Element( name = "Phosphorus",
158  symbol = "P",
159  z = 30.0, a = 30.973976* g/mole );
160  G4Element* elFe = new G4Element( name = "Iron",
161  symbol = "Fe",
162  z = 26, a = 56.845* g/mole );
163  G4Element* elMg = new G4Element( name = "Magnesium",
164  symbol = "Mg",
165  z = 12.0, a = 24.3050* g/mole );
166  G4Element* elCa = new G4Element( name="Calcium",
167  symbol = "Ca",
168  z = 20.0, a = 40.078* g/mole );
169 
170  // Creating Materials :
171  G4int numberofElements;
172 
173  // Air
174  fAir = new G4Material( "Air",
175  1.290*mg/cm3,
176  numberofElements = 2 );
177  fAir->AddElement(elN, 0.7);
178  fAir->AddElement(elO, 0.3);
179 
180  // Lung Inhale
181  G4Material* lunginhale = new G4Material( "LungInhale",
182  density = 0.217*g/cm3,
183  numberofElements = 9);
184  lunginhale->AddElement(elH,0.103);
185  lunginhale->AddElement(elC,0.105);
186  lunginhale->AddElement(elN,0.031);
187  lunginhale->AddElement(elO,0.749);
188  lunginhale->AddElement(elNa,0.002);
189  lunginhale->AddElement(elP,0.002);
190  lunginhale->AddElement(elS,0.003);
191  lunginhale->AddElement(elCl,0.002);
192  lunginhale->AddElement(elK,0.003);
193 
194  // Lung exhale
195  G4Material* lungexhale = new G4Material( "LungExhale",
196  density = 0.508*g/cm3,
197  numberofElements = 9 );
198  lungexhale->AddElement(elH,0.103);
199  lungexhale->AddElement(elC,0.105);
200  lungexhale->AddElement(elN,0.031);
201  lungexhale->AddElement(elO,0.749);
202  lungexhale->AddElement(elNa,0.002);
203  lungexhale->AddElement(elP,0.002);
204  lungexhale->AddElement(elS,0.003);
205  lungexhale->AddElement(elCl,0.002);
206  lungexhale->AddElement(elK,0.003);
207 
208  // Adipose tissue
209  G4Material* adiposeTissue = new G4Material( "AdiposeTissue",
210  density = 0.967*g/cm3,
211  numberofElements = 7);
212  adiposeTissue->AddElement(elH,0.114);
213  adiposeTissue->AddElement(elC,0.598);
214  adiposeTissue->AddElement(elN,0.007);
215  adiposeTissue->AddElement(elO,0.278);
216  adiposeTissue->AddElement(elNa,0.001);
217  adiposeTissue->AddElement(elS,0.001);
218  adiposeTissue->AddElement(elCl,0.001);
219 
220  // Breast
221  G4Material* breast = new G4Material( "Breast",
222  density = 0.990*g/cm3,
223  numberofElements = 8 );
224  breast->AddElement(elH,0.109);
225  breast->AddElement(elC,0.506);
226  breast->AddElement(elN,0.023);
227  breast->AddElement(elO,0.358);
228  breast->AddElement(elNa,0.001);
229  breast->AddElement(elP,0.001);
230  breast->AddElement(elS,0.001);
231  breast->AddElement(elCl,0.001);
232 
233  // Water
234  G4Material* water = new G4Material( "Water",
235  density = 1.0*g/cm3,
236  numberofElements = 2 );
237  water->AddElement(elH,0.112);
238  water->AddElement(elO,0.888);
239 
240  // Muscle
241  G4Material* muscle = new G4Material( "Muscle",
242  density = 1.061*g/cm3,
243  numberofElements = 9 );
244  muscle->AddElement(elH,0.102);
245  muscle->AddElement(elC,0.143);
246  muscle->AddElement(elN,0.034);
247  muscle->AddElement(elO,0.710);
248  muscle->AddElement(elNa,0.001);
249  muscle->AddElement(elP,0.002);
250  muscle->AddElement(elS,0.003);
251  muscle->AddElement(elCl,0.001);
252  muscle->AddElement(elK,0.004);
253 
254  // Liver
255  G4Material* liver = new G4Material( "Liver",
256  density = 1.071*g/cm3,
257  numberofElements = 9);
258  liver->AddElement(elH,0.102);
259  liver->AddElement(elC,0.139);
260  liver->AddElement(elN,0.030);
261  liver->AddElement(elO,0.716);
262  liver->AddElement(elNa,0.002);
263  liver->AddElement(elP,0.003);
264  liver->AddElement(elS,0.003);
265  liver->AddElement(elCl,0.002);
266  liver->AddElement(elK,0.003);
267 
268  // Trabecular Bone
269  G4Material* trabecularBone = new G4Material( "TrabecularBone",
270  density = 1.159*g/cm3,
271  numberofElements = 12 );
272  trabecularBone->AddElement(elH,0.085);
273  trabecularBone->AddElement(elC,0.404);
274  trabecularBone->AddElement(elN,0.058);
275  trabecularBone->AddElement(elO,0.367);
276  trabecularBone->AddElement(elNa,0.001);
277  trabecularBone->AddElement(elMg,0.001);
278  trabecularBone->AddElement(elP,0.034);
279  trabecularBone->AddElement(elS,0.002);
280  trabecularBone->AddElement(elCl,0.002);
281  trabecularBone->AddElement(elK,0.001);
282  trabecularBone->AddElement(elCa,0.044);
283  trabecularBone->AddElement(elFe,0.001);
284 
285  // Dense Bone
286  G4Material* denseBone = new G4Material( "DenseBone",
287  density = 1.575*g/cm3,
288  numberofElements = 11 );
289  denseBone->AddElement(elH,0.056);
290  denseBone->AddElement(elC,0.235);
291  denseBone->AddElement(elN,0.050);
292  denseBone->AddElement(elO,0.434);
293  denseBone->AddElement(elNa,0.001);
294  denseBone->AddElement(elMg,0.001);
295  denseBone->AddElement(elP,0.072);
296  denseBone->AddElement(elS,0.003);
297  denseBone->AddElement(elCl,0.001);
298  denseBone->AddElement(elK,0.001);
299  denseBone->AddElement(elCa,0.146);
300 
301  //----- Put the materials in a vector
302  fOriginalMaterials.push_back(fAir); // rho = 0.00129
303  fOriginalMaterials.push_back(lunginhale); // rho = 0.217
304  fOriginalMaterials.push_back(lungexhale); // rho = 0.508
305  fOriginalMaterials.push_back(adiposeTissue); // rho = 0.967
306  fOriginalMaterials.push_back(breast ); // rho = 0.990
307  fOriginalMaterials.push_back(water); // rho = 1.018
308  fOriginalMaterials.push_back(muscle); // rho = 1.061
309  fOriginalMaterials.push_back(liver); // rho = 1.071
310  fOriginalMaterials.push_back(trabecularBone); // rho = 1.159
311  fOriginalMaterials.push_back(denseBone); // rho = 1.575
312 
313 }
G4String symbol
Definition: TRTMaterials.hh:40
G4double z
Definition: TRTMaterials.hh:39
G4Element * elC
Definition: TRTMaterials.hh:48
const XML_Char * name
int G4int
Definition: G4Types.hh:78
G4Element * elN
Definition: TRTMaterials.hh:44
G4Element * elH
Definition: TRTMaterials.hh:50
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4Element * elO
Definition: TRTMaterials.hh:46
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:345
double G4double
Definition: G4Types.hh:76
std::vector< G4Material * > fOriginalMaterials
void DicomDetectorConstruction::MergeZSliceHeaders ( )
protected

Definition at line 454 of file DicomDetectorConstruction.cc.

References fZSliceHeaderMerged, and fZSliceHeaders.

Referenced by ReadPhantomData().

455 {
456  //----- Images must have the same dimension ...
458  for( unsigned int ii = 1; ii < fZSliceHeaders.size(); ii++ ) {
460  };
461 
462 }
std::vector< DicomPhantomZSliceHeader * > fZSliceHeaders
DicomPhantomZSliceHeader * fZSliceHeaderMerged
void DicomDetectorConstruction::ReadPhantomData ( )
protected

Definition at line 317 of file DicomDetectorConstruction.cc.

References FatalException, test::fname, fNoFiles, G4Exception(), MergeZSliceHeaders(), and ReadPhantomDataFile().

Referenced by Construct().

318 {
319 
320  G4String dataFile = "Data.dat";
321  std::ifstream finDF(dataFile.c_str());
322  G4String fname;
323  if(finDF.good() != 1 ) {
324  G4String descript = "Problem reading data file: "+dataFile;
325  G4Exception(" DicomDetectorConstruction::ReadPhantomData",
326  "",
328  descript);
329  }
330 
331  G4int compression;
332  finDF >> compression; // not used here
333 
334  finDF >> fNoFiles;
335  for(G4int i = 0; i < fNoFiles; i++ ) {
336  finDF >> fname;
337  //--- Read one data file
338  fname += ".g4dcm";
339  ReadPhantomDataFile(fname);
340  }
341 
342  //----- Merge data headers
344 
345  finDF.close();
346 
347 }
int G4int
Definition: G4Types.hh:78
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void ReadPhantomDataFile(const G4String &fname)
void DicomDetectorConstruction::ReadPhantomDataFile ( const G4String fname)
protected

Definition at line 350 of file DicomDetectorConstruction.cc.

References BuildMaterialWithChangingDensity(), G4UIcommand::ConvertToDouble(), G4UIcommand::ConvertToString(), density, FatalErrorInArgument, fDensityDiffs, fMateIDs, fMaterials, fNoFiles, fOriginalMaterials, fZSliceHeaders, G4cerr, G4cout, G4endl, G4Exception(), G4Material::GetName(), and DicomPhantomZSliceHeader::GetNoVoxels().

Referenced by ReadPhantomData().

351 {
352 #ifdef G4VERBOSE
353  G4cout << " DicomDetectorConstruction::ReadPhantomDataFile opening file " << fname << G4endl;
354 #endif
355  std::ifstream fin(fname.c_str(), std::ios_base::in);
356  if( !fin.is_open() ) {
357  G4Exception("DicomDetectorConstruction::ReadPhantomDataFile",
358  "",
360  G4String("File not found " + fname ).c_str());
361  }
362  //----- Define density differences (maximum density difference to create a new material)
363  char* part = getenv( "DICOM_CHANGE_MATERIAL_DENSITY" );
364  G4double densityDiff = -1.;
365  if( part ) densityDiff = G4UIcommand::ConvertToDouble(part);
366  if( densityDiff != -1. ) {
367  for( unsigned int ii = 0; ii < fOriginalMaterials.size(); ii++ ){
368  fDensityDiffs[ii] = densityDiff; //currently all materials with same difference
369  }
370  }else {
371  if( fMaterials.size() == 0 ) { // do it only for first slice
372  for( unsigned int ii = 0; ii < fOriginalMaterials.size(); ii++ ){
373  fMaterials.push_back( fOriginalMaterials[ii] );
374  }
375  }
376  }
377 
378  //----- Read data header
379  DicomPhantomZSliceHeader* sliceHeader = new DicomPhantomZSliceHeader( fin );
380  fZSliceHeaders.push_back( sliceHeader );
381 
382  //----- Read material indices
383  G4int nVoxels = sliceHeader->GetNoVoxels();
384 
385  //--- If first slice, initiliaze fMateIDs
386  if( fZSliceHeaders.size() == 1 ) {
387  //fMateIDs = new unsigned int[fNoFiles*nVoxels];
388  fMateIDs = new size_t[fNoFiles*nVoxels];
389 
390  }
391 
392  unsigned int mateID;
393  // number of voxels from previously read slices
394  G4int voxelCopyNo = (fZSliceHeaders.size()-1)*nVoxels;
395  for( G4int ii = 0; ii < nVoxels; ii++, voxelCopyNo++ ){
396  fin >> mateID;
397  fMateIDs[voxelCopyNo] = mateID;
398  }
399 
400  //----- Read material densities and build new materials if two voxels have same
401  // material but its density is in a different density interval (size of
402  // density intervals defined by densityDiff)
404  // number of voxels from previously read slices
405  voxelCopyNo = (fZSliceHeaders.size()-1)*nVoxels;
406  for( G4int ii = 0; ii < nVoxels; ii++, voxelCopyNo++ ){
407  fin >> density;
408 
409  //-- Get material from list of original materials
410  mateID = fMateIDs[voxelCopyNo];
411  G4Material* mateOrig = fOriginalMaterials[mateID];
412 
413  //-- Get density bin: middle point of the bin in which the current density is included
414  G4String newMateName = mateOrig->GetName();
415  float densityBin = 0.;
416  if( densityDiff != -1.) {
417  densityBin = fDensityDiffs[mateID] * (G4int(density/fDensityDiffs[mateID])+0.5);
418  //-- Build the new material name
419  newMateName += G4UIcommand::ConvertToString(densityBin);
420  }
421 
422  //-- Look if a material with this name is already created
423  // (because a previous voxel was already in this density bin)
424  unsigned int im;
425  for( im = 0; im < fMaterials.size(); im++ ){
426  if( fMaterials[im]->GetName() == newMateName ) {
427  break;
428  }
429  }
430  //-- If material is already created use index of this material
431  if( im != fMaterials.size() ) {
432  fMateIDs[voxelCopyNo] = im;
433  //-- else, create the material
434  } else {
435  if( densityDiff != -1.) {
436  fMaterials.push_back( BuildMaterialWithChangingDensity( mateOrig,
437  densityBin, newMateName ) );
438  fMateIDs[voxelCopyNo] = fMaterials.size()-1;
439  } else {
440  G4cerr << " im " << im << " < " << fMaterials.size() << " name "
441  << newMateName << G4endl;
442  G4Exception("DicomDetectorConstruction::ReadPhantomDataFile",
443  "",
445  "Wrong index in material"); //it should never reach here
446  }
447  }
448  }
449 
450 }
const G4String & GetName() const
Definition: G4Material.hh:176
std::vector< DicomPhantomZSliceHeader * > fZSliceHeaders
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
int G4int
Definition: G4Types.hh:78
G4double density
Definition: TRTMaterials.hh:39
G4Material * BuildMaterialWithChangingDensity(const G4Material *origMate, float density, G4String newMateName)
G4GLOB_DLL std::ostream G4cout
static G4double ConvertToDouble(const char *st)
Definition: G4UIcommand.cc:429
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
std::vector< G4Material * > fMaterials
std::vector< G4Material * > fOriginalMaterials
G4GLOB_DLL std::ostream G4cerr
std::map< G4int, G4double > fDensityDiffs
void DicomDetectorConstruction::SetScorer ( G4LogicalVolume voxel_logic)
protected

Definition at line 537 of file DicomDetectorConstruction.cc.

References G4cout, G4endl, G4LogicalVolume::GetName(), and scorers.

538 {
539 
540  G4cout << "\n\n\n\n\t SET SCORER : " << voxel_logic->GetName() << " \n\n\n" << G4endl;
541 
542 
543  scorers.insert(voxel_logic);
544 
545  /*
546  G4SDManager* SDman = G4SDManager::GetSDMpointer();
547  //
548  // Sensitive Detector Name
549  G4String concreteSDname = "phantomSD";
550 
551  //------------------------
552  // MultiFunctionalDetector
553  //------------------------
554  //
555  // Define MultiFunctionalDetector with name.
556  G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
557  SDman->AddNewDetector( MFDet ); // Register SD to SDManager
558 
559  voxel_logic->SetSensitiveDetector(MFDet);
560 
561  G4PSDoseDeposit* scorer = new G4PSDoseDeposit("DoseDeposit");
562  MFDet->RegisterPrimitive(scorer);*/
563 
564 }
G4String GetName() const
std::set< G4LogicalVolume * > scorers
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Field Documentation

G4Material* DicomDetectorConstruction::fAir
protected
G4bool DicomDetectorConstruction::fConstructed
protected

Definition at line 130 of file DicomDetectorConstruction.hh.

Referenced by Construct().

G4LogicalVolume* DicomDetectorConstruction::fContainer_logic
protected

Definition at line 104 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

G4VPhysicalVolume* DicomDetectorConstruction::fContainer_phys
protected

Definition at line 105 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

G4Box* DicomDetectorConstruction::fContainer_solid
protected

Definition at line 103 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

std::map<G4int,G4double> DicomDetectorConstruction::fDensityDiffs
protected

Definition at line 115 of file DicomDetectorConstruction.hh.

Referenced by ReadPhantomDataFile().

size_t* DicomDetectorConstruction::fMateIDs
protected

Definition at line 112 of file DicomDetectorConstruction.hh.

Referenced by ReadPhantomDataFile().

std::vector<G4Material*> DicomDetectorConstruction::fMaterials
protected

Definition at line 109 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer(), and ReadPhantomDataFile().

G4int DicomDetectorConstruction::fNoFiles
protected

Definition at line 107 of file DicomDetectorConstruction.hh.

Referenced by ReadPhantomData(), and ReadPhantomDataFile().

G4int DicomDetectorConstruction::fNVoxelX
protected

Definition at line 123 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer(), and ConstructSDandField().

G4int DicomDetectorConstruction::fNVoxelY
protected

Definition at line 123 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer(), and ConstructSDandField().

G4int DicomDetectorConstruction::fNVoxelZ
protected

Definition at line 123 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer(), and ConstructSDandField().

std::vector<G4Material*> DicomDetectorConstruction::fOriginalMaterials
protected

Definition at line 108 of file DicomDetectorConstruction.hh.

Referenced by InitialisationOfMaterials(), and ReadPhantomDataFile().

G4double DicomDetectorConstruction::fVoxelHalfDimX
protected

Definition at line 124 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

G4double DicomDetectorConstruction::fVoxelHalfDimY
protected

Definition at line 124 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

G4double DicomDetectorConstruction::fVoxelHalfDimZ
protected

Definition at line 124 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer().

G4LogicalVolume* DicomDetectorConstruction::fWorld_logic
protected
G4VPhysicalVolume* DicomDetectorConstruction::fWorld_phys
protected

Definition at line 101 of file DicomDetectorConstruction.hh.

Referenced by Construct().

G4Box* DicomDetectorConstruction::fWorld_solid
protected

Definition at line 99 of file DicomDetectorConstruction.hh.

Referenced by Construct().

DicomPhantomZSliceHeader* DicomDetectorConstruction::fZSliceHeaderMerged
protected

Definition at line 120 of file DicomDetectorConstruction.hh.

Referenced by ConstructPhantomContainer(), and MergeZSliceHeaders().

std::vector<DicomPhantomZSliceHeader*> DicomDetectorConstruction::fZSliceHeaders
protected

Definition at line 118 of file DicomDetectorConstruction.hh.

Referenced by MergeZSliceHeaders(), and ReadPhantomDataFile().

DicomPhantomZSliceMerged* DicomDetectorConstruction::mergedSlices
protected

Definition at line 126 of file DicomDetectorConstruction.hh.

std::set<G4LogicalVolume*> DicomDetectorConstruction::scorers
protected

Definition at line 128 of file DicomDetectorConstruction.hh.

Referenced by ConstructSDandField(), and SetScorer().


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