Geant4  10.03
DicomPartialDetectorConstruction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: DicomPartialDetectorConstruction.cc 101301 2016-11-14 11:19:22Z gcosmo $
27 //
30 //
31 
32 #include "globals.hh"
33 
35 
36 #include "G4Box.hh"
37 #include "G4LogicalVolume.hh"
38 #include "G4VPhysicalVolume.hh"
39 #include "G4PVPlacement.hh"
40 #include "G4PVParameterised.hh"
41 #include "G4Material.hh"
42 #include "G4Element.hh"
43 #include "G4VisAttributes.hh"
44 #include "G4Colour.hh"
45 #include "G4ios.hh"
47 #include "G4NistManager.hh"
48 #include "G4UIcommand.hh"
49 #include "G4Tubs.hh"
50 
51 #include "G4tgbVolumeMgr.hh"
52 #include "G4tgbMaterialMgr.hh"
53 #include "G4tgrMessenger.hh"
54 #include "G4tgrMessenger.hh"
55 
56 #include "G4SystemOfUnits.hh"
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61  fPartialPhantomParam(0),
62  fNVoxels(0),
63  fDimX(0),
64  fDimY(0),
65  fDimZ(0),
66  fOffsetX(0),
67  fOffsetY(0),
68  fOffsetZ(0),
69  fMateIDs(0)
70 {
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75 {
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 {
82 
83  //----- Build world
84  G4double worldXDimension = 1.*m;
85  G4double worldYDimension = 1.*m;
86  G4double worldZDimension = 1.*m;
87 
88  G4Box* world_box = new G4Box( "WorldSolid",
89  worldXDimension,
90  worldYDimension,
91  worldZDimension );
92 
93  fWorld_logic = new G4LogicalVolume( world_box,
94  fAir,
95  "WorldLogical",
96  0, 0, 0 );
97 
98  G4VPhysicalVolume* world_phys = new G4PVPlacement( 0,
99  G4ThreeVector(0,0,0),
100  "World",
101  fWorld_logic,
102  0,
103  false,
104  0 );
105 
106  ReadPhantomData();
107 
110 
111  return world_phys;
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 {
117  G4String dataFile = "Data.dat";
118  std::ifstream finDF(dataFile.c_str());
119  G4String fname;
120  if(finDF.good() != 1 ) {
121  G4Exception(" DicomPartialDetectorConstruction::ReadPhantomData",
122  "",
124  " Problem reading data file: Data.dat");
125  }
126 
127  G4int compression;
128  finDF >> compression; // not used here
129 
130  G4int numFiles;
131  finDF >> numFiles; // only 1 file supported for the moment
132  if( numFiles != 1 ) {
133  G4Exception("DicomPartialDetectorConstruction::ReadPhantomData",
134  "",
136  "More than 1 DICOM file is not supported");
137  }
138  for(G4int i = 0; i < numFiles; i++ ) {
139  finDF >> fname;
140  //--- Read one data file
141  fname += ".g4pdcm";
142  ReadPhantomDataFile(fname);
143  }
144 
145  finDF.close();
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 const G4String& fname)
151 {
152  // G4String filename = "phantom.g4pdcm";
153 #ifdef G4VERBOSE
154  G4cout << " DicomDetectorConstruction::ReadPhantomDataFile opening file "
155  << fname << G4endl;
156 #endif
157  std::ifstream fin(fname.c_str(), std::ios_base::in);
158  if( !fin.is_open() ) {
159  G4Exception("DicomPartialDetectorConstruction::ReadPhantomDataFile",
160  "",
162  G4String("File not found " + fname).c_str());
163  }
164  G4int nMaterials;
165  fin >> nMaterials;
166  G4String stemp;
167  G4int nmate;
168  for( G4int ii = 0; ii < nMaterials; ii++ ){
169  fin >> nmate >> stemp;
170  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData reading nmate "
171  << ii << " = " << nmate << " mate " << stemp << G4endl;
172  if( ii != nmate )
173  G4Exception("DicomPartialDetectorConstruction::ReadPhantomData",
174  "",
176  "Material number should be in increasing order: \
177  wrong material number ");
178 
179  }
180 
181  fin >> fNVoxelX >> fNVoxelY >> fNVoxelZ;
182  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData nVoxel X/Y/Z "
183  << fNVoxelX << " " << fNVoxelY << " " << fNVoxelZ << G4endl;
184  fin >> fOffsetX >> fDimX;
185  fDimX = (fDimX - fOffsetX)/fNVoxelX;
186  fin >> fOffsetY >> fDimY;
187  fDimY = (fDimY - fOffsetY)/fNVoxelY;
188  fin >> fOffsetZ >> fDimZ;
189  fDimZ = (fDimZ - fOffsetZ)/fNVoxelZ;
190  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData voxelDimX "
191  << fDimX << " fOffsetX " << fOffsetX << G4endl;
192  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData voxelDimY "
193  << fDimY << " fOffsetY " << fOffsetY << G4endl;
194  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData voxelDimZ "
195  << fDimZ << " fOffsetZ " << fOffsetZ << G4endl;
196 
197  //--- Read voxels that are filled
198  fNVoxels = 0;
199  // G4bool* isFilled = new G4bool[fNVoxelX*fNVoxelY*fNVoxelZ];
200  // fFilledIDs = new size_t[fNVoxelZ*fNVoxelY+1];
201  //? fFilledIDs.insert(0);
202  int ifxmin1, ifxmax1;
203  for( G4int iz = 0; iz < fNVoxelZ; iz++ ) {
204  std::map< G4int, G4int > ifmin, ifmax;
205  for( G4int iy = 0; iy < fNVoxelY; iy++ ) {
206  fin >> ifxmin1 >> ifxmax1;
207  // check coherence ...
208 
209  ifmin[iy] = ifxmin1;
210  ifmax[iy] = ifxmax1;
211  G4int ifxdiff = ifxmax1-ifxmin1+1;
212  if( ifxmax1 == -1 && ifxmin1 == -1 ) ifxdiff = 0;
213  fFilledIDs.insert(std::pair<size_t,size_t>(ifxdiff+fNVoxels-1, ifxmin1));
214  G4cout << "DicomPartialDetectorConstruction::ReadPhantomData insert "
215  << " FilledIDs " << ifxdiff+fNVoxels-1 << " min " << ifxmin1
216  << " N= " << fFilledIDs.size() << G4endl;
217  //filledIDs[iz*fNVoxelY+iy+1] = ifxmax1-ifxmin1+1;
218  for( G4int ix = 0; ix < fNVoxelX; ix++ ) {
219  if( ix >= G4int(ifxmin1) && ix <= G4int(ifxmax1) ) {
220  fNVoxels++;
221  } else {
222  }
223  }
224  }
225  fFilledMins[iz] = ifmin;
226  fFilledMaxs[iz] = ifmax;
227  }
228 
229  //--- Read material IDs
230  G4int mateID1;
231  fMateIDs = new size_t[fNVoxelX*fNVoxelY*fNVoxelZ];
232  G4int copyNo = 0;
233  for( G4int iz = 0; iz < fNVoxelZ; iz++ ) {
234  std::map< G4int, G4int > ifmin = fFilledMins[iz];
235  std::map< G4int, G4int > ifmax = fFilledMaxs[iz];
236  for( G4int iy = 0; iy < fNVoxelY; iy++ ) {
237  ifxmin1 = ifmin[iy];
238  ifxmax1 = ifmax[iy];
239  for( G4int ix = 0; ix < fNVoxelX; ix++ ) {
240  if( ix >= G4int(ifxmin1) && ix <= G4int(ifxmax1) ) {
241  fin >> mateID1;
242  if( mateID1 < 0 || mateID1 >= nMaterials ) {
243  G4Exception("DicomPartialDetectorConstruction::ReadPhantomData",
244  "",
246  G4String("Wrong index in phantom file: It should be between 0 and "
247  + G4UIcommand::ConvertToString(nMaterials-1)
248  + ", while it is "
249  + G4UIcommand::ConvertToString(mateID1)).c_str());
250  }
251  fMateIDs[copyNo] = mateID1;
252  copyNo++;
253  }
254  }
255  }
256  }
257 
259 
260  fin.close();
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265 {
266  std::map<G4int, std::pair<G4double,G4double> > densiMinMax;
267  std::map<G4int, std::pair<G4double,G4double> >::iterator mpite;
268  for( size_t ii = 0; ii < fOriginalMaterials.size(); ii++ ){
269  densiMinMax[ii] = std::pair<G4double,G4double>(DBL_MAX,-DBL_MAX);
270  }
271 
272  //----- Define density differences (maximum density difference to create
273  // a new material)
274  char* part = getenv( "DICOM_CHANGE_MATERIAL_DENSITY" );
275  G4double densityDiff = -1.;
276  if( part ) densityDiff = G4UIcommand::ConvertToDouble(part);
277  std::map<G4int,G4double> densityDiffs;
278  if( densityDiff != -1. ) {
279  for( size_t ii = 0; ii < fOriginalMaterials.size(); ii++ ){
280  densityDiffs[ii] = densityDiff; //currently all materials
281  //with same step
282  }
283  }else {
284  if( fMaterials.size() == 0 ) { // do it only for first slice
285  for( unsigned int ii = 0; ii < fOriginalMaterials.size(); ii++ ){
286  fMaterials.push_back( fOriginalMaterials[ii] );
287  }
288  }
289  }
290  // densityDiffs[0] = 0.0001; //fAir
291 
292  //--- Calculate the average material density for each material/density bin
293  std::map< std::pair<G4Material*,G4int>, matInfo* > newMateDens;
294 
295  G4double dens1;
296  size_t ifxmin1, ifxmax1;
297 
298  //---- Read the material densities
299  G4int copyNo = 0;
300  for( G4int iz = 0; iz < fNVoxelZ; iz++ ) {
301  std::map< G4int, G4int > ifmin = fFilledMins[iz];
302  std::map< G4int, G4int > ifmax = fFilledMaxs[iz];
303  for( G4int iy = 0; iy < fNVoxelY; iy++ ) {
304  ifxmin1 = ifmin[iy];
305  ifxmax1 = ifmax[iy];
306  for( G4int ix = 0; ix < fNVoxelX; ix++ ) {
307  if( ix >= G4int(ifxmin1) && ix <= G4int(ifxmax1) ) {
308  fin >> dens1;
309  //--- store the minimum and maximum density for each material
310  //(just for printing)
311  mpite = densiMinMax.find( fMateIDs[copyNo] );
312  if( dens1 < (*mpite).second.first ) (*mpite).second.first = dens1;
313  if( dens1 > (*mpite).second.second ) (*mpite).second.second = dens1;
314 
315  //--- Get material from original list of material in file
316  int mateID = fMateIDs[copyNo];
317  G4Material* mate = fOriginalMaterials[mateID];
318  // G4cout << copyNo << " mateID " << mateID << G4endl;
319  //--- Check if density is equal to the original material density
320  if( std::fabs(dens1 - mate->GetDensity()/g*cm3 ) < 1.e-9 ) continue;
321 
322  //--- Build material name with fOriginalMaterials name + density
323  // float densityBin = densityDiffs[mateID] *
324  // (G4int(dens1/densityDiffs[mateID])+0.5);
325  G4String newMateName = mate->GetName();
326 
327  G4int densityBin = 0;
328  // G4cout << " densityBin " << densityBin << " "
329  // << dens1 << " " <<densityDiffs[mateID] << G4endl;
330 
331  if( densityDiff != -1.) {
332  densityBin = (G4int(dens1/densityDiffs[mateID]));
333  newMateName =
334  mate->GetName()+G4UIcommand::ConvertToString(densityBin);
335  }
336 
337  //--- Look if it is the first voxel with this material/densityBin
338  std::pair<G4Material*,G4int> matdens(mate, densityBin );
339 
340  std::map< std::pair<G4Material*,G4int>, matInfo* >::iterator mppite
341  = newMateDens.find( matdens );
342  if( mppite != newMateDens.end() ){
343  matInfo* mi = (*mppite).second;
344  mi->fSumdens += dens1;
345  mi->fNvoxels++;
346  fMateIDs[copyNo] = fOriginalMaterials.size()-1 + mi->fId;
347  // G4cout << copyNo << " mat new again "
348  //<< fOriginalMaterials.size()-1 + mi->id << " " << mi->id << G4endl;
349  } else {
350  matInfo* mi = new matInfo;
351  mi->fSumdens = dens1;
352  mi->fNvoxels = 1;
353  mi->fId = newMateDens.size()+1;
354  newMateDens[matdens] = mi;
355  fMateIDs[copyNo] = fOriginalMaterials.size()-1 + mi->fId;
356  // G4cout << copyNo << " mat new first "
357  // << fOriginalMaterials.size()-1 + mi->fId << G4endl;
358  }
359  copyNo++;
360  // G4cout << ix << " " << iy << " " << iz
361  // << " filling fMateIDs " << copyNo << " = " << atoi(cid)-1 << G4endl;
362  // fMateIDs[copyNo] = atoi(cid)-1;
363  }
364  }
365  }
366  }
367 
368  //----- Build the list of phantom materials that go to Parameterisation
369  //--- Add original materials
370  std::vector<G4Material*>::const_iterator mimite;
371  for( mimite = fOriginalMaterials.begin(); mimite != fOriginalMaterials.end();
372  mimite++ ){
373  fPhantomMaterials.push_back( (*mimite) );
374  }
375  //
376  //---- Build and add new materials
377  std::map< std::pair<G4Material*,G4int>, matInfo* >::iterator mppite;
378  for( mppite= newMateDens.begin(); mppite != newMateDens.end(); mppite++ ){
379  G4double averdens = (*mppite).second->fSumdens/(*mppite).second->fNvoxels;
380  G4double saverdens = G4int(1000.001*averdens)/1000.;
381  G4cout << "GmReadPhantomGeometry::ReadVoxelDensities AVER DENS "
382  << averdens << " -> " << saverdens << " -> "
383  << G4int(1000*averdens) << " " << G4int(1000*averdens)/1000 << " "
384  << G4int(1000*averdens)/1000. << G4endl;
385 
386  G4cout << "GmReadPhantomGeometry::ReadVoxelDensities AVER DENS "
387  << averdens << " -> " << saverdens << " -> "
388  << G4UIcommand::ConvertToString(saverdens) << " Nvoxels "
389  <<(*mppite).second->fNvoxels << G4endl;
390  G4String mateName = ((*mppite).first).first->GetName() + "_" +
391  G4UIcommand::ConvertToString(saverdens);
393  (*mppite).first.first,
394  averdens, mateName ) );
395  }
396 
397 }
398 
399 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401 {
402  //build a clinder that encompass the partial phantom built in the example
403  // /dicom/intersectWithUserVolume 0. 0. 0. 0.*deg 0. 0. TUBE 0. 200. 100.
404  //---- Extract number of voxels and voxel dimensions
405 
406  G4String contname= "PhantomContainer";
407 
408  //----- Define the volume that contains all the voxels
409  G4Tubs* container_tube = new G4Tubs(contname,0., 200., 100., 0., 360*deg);
410 
412  new G4LogicalVolume( container_tube,
413  fAir,
414  contname,
415  0, 0, 0 );
416 
417  G4ThreeVector posCentreVoxels(0.,0.,0.);
418  //G4cout << " PhantomContainer posCentreVoxels " << posCentreVoxels << G4endl;
420 
422  new G4PVPlacement(rotm, // rotation
423  posCentreVoxels,
424  fContainer_logic, // The logic volume
425  contname, // Name
426  fWorld_logic, // Mother
427  false, // No op. bool.
428  1); // Copy number
429 
430 }
431 
432 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
434 {
435  G4String OptimAxis = "kXAxis";
436  G4bool bSkipEqualMaterials = 0;
437  G4int regStructureID = 2;
438 
439  G4ThreeVector posCentreVoxels(0.,0.,0.);
440 
441  //----- Build phantom
442  G4String voxelName = "phantom";
443  G4Box* phantom_solid = new G4Box(voxelName, fDimX/2., fDimY/2., fDimZ/2. );
444  G4LogicalVolume* phantom_logic =
445  new G4LogicalVolume( phantom_solid,
446  fAir,
447  voxelName,
448  0, 0, 0 );
449  G4bool pVis = 0;
450  if( !pVis ) {
451  G4VisAttributes* visAtt = new G4VisAttributes( FALSE );
452  phantom_logic->SetVisAttributes( visAtt );
453  }
454 
455  G4double theSmartless = 0.5;
456  // fContainer_logic->SetSmartless( theSmartless );
457  phantom_logic->SetSmartless( theSmartless );
458 
460 
465 
467 
469 
471 
472  // G4cout << " Number of Materials " << fPhantomMaterials.size() << G4endl;
473  // G4cout << " SetMaterialIndices(0) " << fMateIDs[0] << G4endl;
474 
475  G4PVParameterised * phantom_phys = 0;
476  if( OptimAxis == "kUndefined" ) {
477  phantom_phys =
478  new G4PVParameterised(voxelName,phantom_logic,fContainer_logic,
480  } else if( OptimAxis == "kXAxis" ) {
481  // G4cout << " optim kX " << G4endl;
482  phantom_phys =
483  new G4PVParameterised(voxelName,phantom_logic,fContainer_logic,
484  // kXAxis, fNVoxelX*fNVoxelY*fNVoxelZ, fPartialPhantomParam);
486  fPartialPhantomParam->SetSkipEqualMaterials(bSkipEqualMaterials);
487  } else {
488  G4Exception("GmReadPhantomGeometry::ConstructPhantom",
489  "",
491  G4String("Wrong argument to parameter \
492  GmReadPhantomGeometry:Phantom:OptimAxis: \
493  Only allowed 'kUndefined' or 'kXAxis', it is: "+OptimAxis).c_str());
494  }
495 
496  phantom_phys->SetRegularStructureId(regStructureID);
497 
498 }
std::map< G4int, std::map< G4int, G4int > > fFilledMins
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
Definition: G4Box.hh:64
const G4String & GetName() const
Definition: G4Material.hh:178
Definition: G4Tubs.hh:85
G4double GetDensity() const
Definition: G4Material.hh:180
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
static constexpr double g
Definition: G4SIunits.hh:183
int G4int
Definition: G4Types.hh:78
void SetFilledIDs(std::multimap< G4int, G4int > fid)
G4PartialPhantomParameterisation * fPartialPhantomParam
void ReadVoxelDensitiesPartial(std::ifstream &fin, std::map< G4int, std::map< G4int, G4int > > ifxmin, std::map< G4int, std::map< G4int, G4int > > ifxmax)
G4Material * BuildMaterialWithChangingDensity(const G4Material *origMate, float density, G4String newMateName)
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
void SetVoxelDimensions(G4double halfx, G4double halfy, G4double halfz)
Definition of the DicomPartialDetectorConstruction class.
bool G4bool
Definition: G4Types.hh:79
void SetFilledMins(std::map< G4int, std::map< G4int, G4int > > fmins)
#define FALSE
Definition: globals.hh:52
Dicom detector construction.
static G4double ConvertToDouble(const char *st)
Definition: G4UIcommand.cc:443
static constexpr double cm3
Definition: G4SIunits.hh:121
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void SetRegularStructureId(G4int Code)
void SetSmartless(G4double s)
std::map< G4int, std::map< G4int, G4int > > fFilledMaxs
void SetNoVoxel(size_t nx, size_t ny, size_t nz)
#define G4endl
Definition: G4ios.hh:61
void SetMaterialIndices(size_t *matInd)
double G4double
Definition: G4Types.hh:76
void SetMaterials(std::vector< G4Material * > &mates)
std::vector< G4Material * > fMaterials
static constexpr double deg
Definition: G4SIunits.hh:152
std::vector< G4Material * > fOriginalMaterials
#define DBL_MAX
Definition: templates.hh:83
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSkipEqualMaterials(G4bool skip)