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