Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DicomIntersectVolume.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 //
28 
29 #include "DicomIntersectVolume.hh"
30 
31 #include "G4UIcmdWithAString.hh"
32 #include "G4LogicalVolume.hh"
33 #include "G4LogicalVolumeStore.hh"
34 #include "G4VPhysicalVolume.hh"
35 #include "G4PhysicalVolumeStore.hh"
36 #include "G4VSolid.hh"
37 #include "G4tgrSolid.hh"
38 #include "G4tgbVolume.hh"
39 #include "G4Material.hh"
41 #include "G4PVParameterised.hh"
42 #include "G4UIcommand.hh"
43 #include "G4SystemOfUnits.hh"
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47  : G4UImessenger(),
48  fSolid(0),
49  fPhantomVolume(0),
50  fRegularParam(0),
51  fPhantomMinusCorner(),
52  fVoxelIsInside(0)
53 {
54  fUserVolumeCmd = new G4UIcmdWithAString("/dicom/intersectWithUserVolume",this);
55  fUserVolumeCmd->SetGuidance("Intersects a phantom with a user-defined volume and outputs the voxels that are totally inside the intersection as a new phantom file. It must have the parameters: POS_X POS_Y POS_Z ANG_X ANG_Y ANG_Z SOLID_TYPE SOLID_PARAM_1 (SOLID_PARAM_2 ...)");
56  fUserVolumeCmd->SetParameterName("choice",true);
57  fUserVolumeCmd->AvailableForStates(G4State_Idle);
58 
59  fG4VolumeCmd = new G4UIcmdWithAString("/dicom/intersectWithUserVolume",this);
60  fG4VolumeCmd->SetGuidance("Intersects a phantom with a user-defined volume and outputs the voxels that are totally inside the intersection as a new phantom file. It must have the parameters: POS_X POS_Y POS_Z ANG_X ANG_Y ANG_Z SOLID_TYPE SOLID_PARAM_1 (SOLID_PARAM_2 ...)");
61  fG4VolumeCmd->SetParameterName("choice",true);
62  fG4VolumeCmd->AvailableForStates(G4State_Idle);
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 {
68  if (fUserVolumeCmd) delete fUserVolumeCmd;
69  if (fG4VolumeCmd) delete fG4VolumeCmd;
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74  G4String newValues)
75 {
76  G4AffineTransform theVolumeTransform;
77 
78  if (command == fUserVolumeCmd) {
79 
80  std::vector<G4String> params = GetWordsInString( newValues );
81  if( params.size() < 8 ) {
82  G4Exception("DicomIntersectVolume::SetNewValue",
83  " There must be at least 8 parameter: SOLID_TYPE POS_X POS_Y POS_Z ANG_X ANG_Y ANG_Z SOLID_PARAM_1 (SOLID_PARAM_2 ...)",
85  G4String("Number of parameters given = " + G4UIcommand::ConvertToString( G4int(params.size()) )).c_str());
86 
87  }
88 
89  //----- Build G4VSolid
90  BuildUserSolid(params);
91 
92  //----- Calculate volume inverse 3D transform
94  G4RotationMatrix* rotmat = new G4RotationMatrix;
95  std::vector<double> angles;
96  rotmat->rotateX( G4UIcommand::ConvertToDouble(params[3]) );
97  rotmat->rotateY( G4UIcommand::ConvertToDouble(params[4]) );
98  rotmat->rotateY( G4UIcommand::ConvertToDouble(params[5]) );
99  theVolumeTransform = G4AffineTransform( rotmat, pos ).Invert();
100 
101  } else if (command == fG4VolumeCmd) {
102  std::vector<G4String> params = GetWordsInString( newValues );
103  if( params.size() !=1 ) G4Exception("DicomIntersectVolume::SetNewValue",
104  "",
106  G4String("Command: "+ command->GetCommandPath() + "/" + command->GetCommandName() + " " + newValues + " needs 1 argument: VOLUME_NAME").c_str());
107 
108  //----- Build G4VSolid
109  BuildG4Solid(params);
110 
111  //----- Calculate volume inverse 3D transform
112  G4VPhysicalVolume* pv = GetPhysicalVolumes( params[0], 1, 1)[0];
113 
114  theVolumeTransform = G4AffineTransform( pv->GetFrameRotation(), pv->GetFrameTranslation() );
115  }
116 
117  //----- Calculate relative phantom - volume 3D transform
118  G4PhantomParameterisation* thePhantomParam = GetPhantomParam(true);
119 
120  G4RotationMatrix* rotph = new G4RotationMatrix(); // assumes the phantom mother is not rotated !!!
121  G4AffineTransform thePhantomTransform( rotph, G4ThreeVector() ); // assumes the phantom mother is not translated !!!
122 
123  G4AffineTransform theTransform = theVolumeTransform*thePhantomTransform;
124 
125  G4ThreeVector axisX( 1., 0., 0. );
126  G4ThreeVector axisY( 0., 1., 0. );
127  G4ThreeVector axisZ( 0., 0., 1. );
128  theTransform.ApplyAxisTransform(axisX);
129  theTransform.ApplyAxisTransform(axisY);
130  theTransform.ApplyAxisTransform(axisZ);
131 
132  //----- Write phantom header
133  G4String thePhantomFileName = "phantom.g4pdcm";
134  fout.open(thePhantomFileName);
135  std::vector<G4Material*> materials = thePhantomParam->GetMaterials();
136  fout << materials.size() << G4endl;
137  for( unsigned int ii = 0; ii < materials.size(); ii++ ) {
138  fout << ii << " " << materials[ii]->GetName() << G4endl;
139  }
140 
141  //----- Loop to pantom voxels
142  G4int nx = thePhantomParam->GetNoVoxelX();
143  G4int ny = thePhantomParam->GetNoVoxelY();
144  G4int nz = thePhantomParam->GetNoVoxelZ();
145  G4int nxy = nx*ny;
146  fVoxelIsInside = new G4bool[nx*ny*nz];
147  G4double voxelHalfWidthX = thePhantomParam->GetVoxelHalfX();
148  G4double voxelHalfWidthY = thePhantomParam->GetVoxelHalfY();
149  G4double voxelHalfWidthZ = thePhantomParam->GetVoxelHalfZ();
150 
151  //----- Write phantom dimensions and limits
152  fout << nx << " " << ny << " " << nz << G4endl;
153  fout << -voxelHalfWidthX*nx+thePhantomTransform.NetTranslation().x() << " "
154  << voxelHalfWidthX*nx+thePhantomTransform.NetTranslation().x() << G4endl;
155  fout << -voxelHalfWidthY*ny+thePhantomTransform.NetTranslation().y() << " "
156  << voxelHalfWidthY*ny+thePhantomTransform.NetTranslation().y() << G4endl;
157  fout << -voxelHalfWidthZ*nz+thePhantomTransform.NetTranslation().z() << " "
158  << voxelHalfWidthZ*nz+thePhantomTransform.NetTranslation().z() << G4endl;
159 
160  for( G4int iz = 0; iz < nz; iz++ ){
161  for( G4int iy = 0; iy < ny; iy++) {
162 
163  G4bool bPrevVoxelInside = true;
164  G4bool b1VoxelFoundInside = false;
165  G4int firstVoxel = -1;
166  G4int lastVoxel = -1;
167  for(G4int ix = 0; ix < nx; ix++ ){
168  G4ThreeVector voxelCentre( (-nx+ix*2+1)*voxelHalfWidthX, (-ny+iy*2+1)*voxelHalfWidthY, (-nz+iz*2+1)*voxelHalfWidthZ);
169  theTransform.ApplyPointTransform(voxelCentre);
170  G4bool bVoxelIsInside = true;
171  for( G4int ivx = -1; ivx <= 1; ivx+=2 ) {
172  for( G4int ivy = -1; ivy <= 1; ivy+=2 ){
173  for( G4int ivz = -1; ivz <= 1; ivz+=2 ) {
174  G4ThreeVector voxelPoint = voxelCentre + ivx*voxelHalfWidthX*axisX + ivy*voxelHalfWidthY*axisY + ivz*voxelHalfWidthZ*axisZ;
175  if( fSolid->Inside( voxelPoint ) == kOutside ) {
176  bVoxelIsInside = false;
177  break;
178  } else {
179  }
180  }
181  if( !bVoxelIsInside ) break;
182  }
183  if( !bVoxelIsInside ) break;
184  }
185 
186  G4int copyNo = ix + nx*iy + nxy*iz;
187  if( bVoxelIsInside ) {
188  if( !bPrevVoxelInside ) {
189  G4Exception("DicomIntersectVolume::SetNewValue",
190  "",
192  "Volume cannot intersect phantom in discontiguous voxels, please use other voxel");
193  }
194  if( !b1VoxelFoundInside ) {
195  firstVoxel = ix;
196  b1VoxelFoundInside = true;
197  }
198  lastVoxel = ix;
199  fVoxelIsInside[copyNo] = true;
200  } else {
201  fVoxelIsInside[copyNo] = false;
202  }
203 
204  }
205  fout << firstVoxel << " " << lastVoxel << G4endl;
206  }
207  }
208 
209  //----- Now write the materials
210  for( G4int iz = 0; iz < nz; iz++ ){
211  for( G4int iy = 0; iy < ny; iy++) {
212  G4bool b1xFound = false;
213  for(G4int ix = 0; ix < nx; ix++ ){
214  size_t copyNo = ix + ny*iy + nxy*iz;
215  // fout << " iz " << iz << " i " << iy << " ix " << ix << G4endl;
216  if( fVoxelIsInside[copyNo] ) {
217  fout << thePhantomParam->GetMaterialIndex(copyNo)<< " ";
218  b1xFound = true;
219  }
220  }
221  if(b1xFound ) fout << G4endl;
222  }
223  }
224 
225  // Now write densities
226  for( G4int iz = 0; iz < nz; iz++ ){
227  for( G4int iy = 0; iy < ny; iy++) {
228  G4bool b1xFound = false;
229  for(G4int ix = 0; ix < nx; ix++ ){
230  size_t copyNo = ix + ny*iy + nxy*iz;
231  if( fVoxelIsInside[copyNo] ) {
232  fout << thePhantomParam->GetMaterial(copyNo)->GetDensity()/g*cm3 << " ";
233  b1xFound = true;
234  }
235  }
236  if(b1xFound ) fout << G4endl;
237  }
238  }
239 
240 }
241 
242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243 void DicomIntersectVolume::BuildUserSolid( std::vector<G4String> params )
244 {
245  for( G4int ii = 0; ii < 6; ii++ ) params.erase( params.begin() ); // take otu position and rotation angles
246  params.insert( params.begin(), ":SOLID");
247  params.insert( params.begin(), params[1] );
248  G4tgrSolid* tgrSolid = new G4tgrSolid(params);
249  G4tgbVolume* tgbVolume = new G4tgbVolume();
250  fSolid = tgbVolume->FindOrConstructG4Solid( tgrSolid );
251 
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255 void DicomIntersectVolume::BuildG4Solid( std::vector<G4String> params )
256 {
257  fSolid = GetLogicalVolumes( params[0], 1, 1)[0]->GetSolid();
258 
259 }
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262 G4PhantomParameterisation* DicomIntersectVolume::GetPhantomParam(G4bool bMustExist)
263 {
264  G4PhantomParameterisation* paramreg = 0;
265 
267  std::vector<G4VPhysicalVolume*>::iterator cite;
268  for( cite = pvs->begin(); cite != pvs->end(); cite++ ) {
269  // G4cout << " PV " << (*cite)->GetName() << " " << (*cite)->GetTranslation() << G4endl;
270  if( IsPhantomVolume( *cite ) ) {
271  const G4PVParameterised* pvparam = static_cast<const G4PVParameterised*>(*cite);
272  G4VPVParameterisation* param = pvparam->GetParameterisation();
273  // if( static_cast<const G4PhantomParameterisation*>(param) ){
274  // if( static_cast<const G4PhantomParameterisation*>(param) ){
275  // G4cout << "G4PhantomParameterisation volume found " << (*cite)->GetName() << G4endl;
276  paramreg = static_cast<G4PhantomParameterisation*>(param);
277  }
278  }
279 
280  if( !paramreg && bMustExist )
281  G4Exception("DicomIntersectVolume::GetPhantomParam",
282  "",
284  " No G4PhantomParameterisation found ");
285 
286  return paramreg;
287 
288 }
289 
290 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
291 std::vector<G4VPhysicalVolume*> DicomIntersectVolume::GetPhysicalVolumes( const G4String& name, bool exists, G4int nVols )
292 {
293  std::vector<G4VPhysicalVolume*> vvolu;
294  std::string::size_type ial = name.rfind(":");
295  G4String volname = "";
296  G4int volcopy = 0;
297  if( ial != std::string::npos ) {
298  std::string::size_type ial2 = name.rfind(":",ial-1);
299  if( ial2 != std::string::npos ) {
300  G4Exception("DicomIntersectVolume::GetPhysicalVolumes",
301  "",
303  G4String("Name corresponds to a touchable " + name).c_str());
304  }else {
305  volname = name.substr( 0, ial );
306  volcopy = G4UIcommand::ConvertToInt( name.substr( ial+1, name.length() ).c_str() );
307  }
308  } else {
309  volname = name;
310  volcopy = -1;
311  }
312 
314  std::vector<G4VPhysicalVolume*>::iterator citepv;
315  for( citepv = pvs->begin(); citepv != pvs->end(); citepv++ ) {
316  if( volname == (*citepv)->GetName()
317  && ( (*citepv)->GetCopyNo() == volcopy || -1 == volcopy ) ){
318  vvolu.push_back( *citepv );
319  }
320  }
321 
322  //----- Check that at least one volume was found
323  if( vvolu.size() == 0 ) {
324  if(exists) {
325  G4Exception(" DicomIntersectVolume::GetPhysicalVolumes",
326  "",
328  G4String("No physical volume found with name " + name).c_str());
329  } else {
330  G4cerr << "!!WARNING: DicomIntersectVolume::GetPhysicalVolumes: no physical volume found with name " << name << G4endl;
331  }
332  }
333 
334  if( nVols != -1 && G4int(vvolu.size()) != nVols ) {
335  G4Exception("DicomIntersectVolume::GetLogicalVolumes:",
336  "Wrong number of physical volumes found",
338  ("Number of physical volumes " + G4UIcommand::ConvertToString(G4int(vvolu.size())) + ", requesting " + G4UIcommand::ConvertToString(nVols)).c_str());
339  }
340 
341  return vvolu;
342 }
343 
344 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345 G4bool DicomIntersectVolume::IsPhantomVolume( G4VPhysicalVolume* pv )
346 {
347  EAxis axis;
348  G4int nReplicas;
349  G4double width,offset;
350  G4bool consuming;
351  pv->GetReplicationData(axis,nReplicas,width,offset,consuming);
352  EVolume type = (consuming) ? kReplica : kParameterised;
353  if( type == kParameterised && pv->GetRegularStructureId() == 1 ) {
354  return TRUE;
355  } else {
356  return FALSE;
357  }
358 
359 }
360 
361 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
362 std::vector<G4LogicalVolume*> DicomIntersectVolume::GetLogicalVolumes( const G4String& name, bool exists, G4int nVols )
363 {
364  // G4cout << "GetLogicalVolumes " << name << " " << exists << G4endl;
365  std::vector<G4LogicalVolume*> vvolu;
366  G4int ial = name.rfind(":");
367  if( ial != -1 ) {
368  G4Exception("DicomIntersectVolume::GetLogicalVolumes",
369  "",
371  G4String("Name corresponds to a touchable or physcal volume" + name).c_str());
372  }
373 
375  std::vector<G4LogicalVolume*>::iterator citelv;
376  for( citelv = lvs->begin(); citelv != lvs->end(); citelv++ ) {
377  if( name == (*citelv)->GetName() ) {
378  vvolu.push_back( *citelv );
379  }
380  }
381 
382  //----- Check that at least one volume was found
383  if( vvolu.size() == 0 ) {
384  if(exists) {
385  G4Exception("DicomIntersectVolume::GetLogicalVolumes:","",FatalErrorInArgument,("no logical volume found with name " + name).c_str());
386  } else {
387  G4Exception("DicomIntersectVolume::GetLogicalVolumes:","",JustWarning,("no logical volume found with name " + name).c_str());
388  }
389  }
390 
391  if( nVols != -1 && G4int(vvolu.size()) != nVols ) {
392  G4Exception("DicomIntersectVolume::GetLogicalVolumes:",
393  "Wrong number of logical volumes found",
395  ("Number of logical volumes " + G4UIcommand::ConvertToString(G4int(vvolu.size())) + ", requesting " + G4UIcommand::ConvertToString(nVols)).c_str());
396  }
397 
398  return vvolu;
399 
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403 std::vector<G4String> DicomIntersectVolume::GetWordsInString( const G4String& stemp)
404 {
405  std::vector<G4String> wordlist;
406 
407  //---------- Read a line of file:
408  //----- Clear wordlist
409  G4int ii;
410  const char* cstr = stemp.c_str();
411  int siz = stemp.length();
412  int istart = 0;
413  int nQuotes = 0;
414  bool lastIsBlank = false;
415  bool lastIsQuote = false;
416  for( ii = 0; ii < siz; ii++ ){
417  if(cstr[ii] == '\"' ){
418  if( lastIsQuote ){
419  G4Exception("GmGenUtils:GetWordsFromString","",FatalException, ("There cannot be two quotes together " + stemp).c_str() );
420  }
421  if( nQuotes%2 == 1 ){
422  //close word
423  wordlist.push_back( stemp.substr(istart,ii-istart) );
424  // G4cout << "GetWordsInString new word0 " << wordlist[wordlist.size()-1] << " istart " << istart << " ii " << ii << G4endl;
425  } else {
426  istart = ii+1;
427  }
428  nQuotes++;
429  lastIsQuote = true;
430  lastIsBlank = false;
431  } else if(cstr[ii] == ' ' ){
432  // G4cout << "GetWordsInString blank nQuotes " << nQuotes << " lastIsBlank " << lastIsBlank << G4endl;
433  if( nQuotes%2 == 0 ){
434  if( !lastIsBlank && !lastIsQuote ) {
435  wordlist.push_back( stemp.substr(istart,ii-istart) );
436  // G4cout << "GetWordsInString new word1 " << wordlist[wordlist.size()-1] << " lastIsBlank " << lastIsBlank << G4endl;
437  }
438 
439  istart = ii+1;
440  lastIsQuote = false;
441  lastIsBlank = true;
442  }
443  } else {
444  if( ii == siz-1 ) {
445  wordlist.push_back( stemp.substr(istart,ii-istart+1) );
446  // G4cout << "GetWordsInString new word2 " << wordlist[wordlist.size()-1] << " istart " << istart << G4endl;
447  }
448  lastIsQuote = false;
449  lastIsBlank = false;
450  }
451  }
452  if( nQuotes%2 == 1 ) {
453  G4Exception("GmGenUtils:GetWordsFromString","",FatalException, ("unbalanced quotes in line " + stemp).c_str() );
454  }
455 
456  return wordlist;
457 }
458