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