www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
MeshBaseImageSampler Class Reference

A helper class for reading and sampling images using VTK, but with a Meshbase object. More...

#include <MeshBaseImageSampler.h>

Inheritance diagram for MeshBaseImageSampler:
[legend]

Public Member Functions

 MeshBaseImageSampler (const InputParameters &parameters)
 Constructor. More...
 
virtual Real sample (const Point &p)
 Return the pixel value for the given point. More...
 
virtual void setupImageSampler (MeshBase &mesh)
 Perform initialization of image data. More...
 
std::string fileSuffix ()
 
const std::vector< std::string > & filenames ()
 

Protected Member Functions

void vtkShiftAndScale ()
 Apply image re-scaling using the vtkImageShiftAndRescale object. More...
 
void vtkThreshold ()
 Perform thresholding. More...
 
void vtkMagnitude ()
 Convert the image to greyscale. More...
 
void vtkFlip ()
 Perform image flipping. More...
 
void errorCheck ()
 

Protected Attributes

int _status
 
std::string _file_suffix
 
std::vector< std::string > _filenames
 

Private Member Functions

vtkSmartPointer< vtkImageFlip > imageFlip (const int &axis)
 Helper method for flipping image. More...
 

Private Attributes

vtkSmartPointer< vtkStringArray > _files
 List of file names to extract data. More...
 
vtkImageData * _data
 Complete image data. More...
 
vtkAlgorithmOutput * _algorithm
 VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers... More...
 
vtkSmartPointer< vtkImageReader2 > _image
 Complete image data. More...
 
vtkSmartPointer< vtkImageThreshold > _image_threshold
 Pointer to thresholding filter. More...
 
vtkSmartPointer< vtkImageShiftScale > _shift_scale_filter
 Pointer to the shift and scaling filter. More...
 
vtkSmartPointer< vtkImageMagnitude > _magnitude_filter
 Pointer to the magnitude filter. More...
 
vtkSmartPointer< vtkImageFlip > _flip_filter
 Pointers to image flipping filter. May be used for x, y, or z. More...
 
Point _origin
 Origin of image. More...
 
std::vector< int > _dims
 Pixel dimension of image. More...
 
Point _physical_dims
 Physical dimensions of image. More...
 
std::vector< double > _voxel
 Physical pixel size. More...
 
unsigned int _component
 Component to extract. More...
 
BoundingBox _bounding_box
 Bounding box for testing points. More...
 
const InputParameters_is_pars
 Parameters for interface. More...
 
ConsoleStream _is_console
 Create a console stream object for this helper class. More...
 

Detailed Description

A helper class for reading and sampling images using VTK, but with a Meshbase object.

Definition at line 52 of file MeshBaseImageSampler.h.

Constructor & Destructor Documentation

◆ MeshBaseImageSampler()

MeshBaseImageSampler::MeshBaseImageSampler ( const InputParameters parameters)

Constructor.

Use this object as an interface, being sure to also add the parameters to the child class.

See also
ImageFunction

Definition at line 57 of file MeshBaseImageSampler.C.

58  : FileRangeBuilder(parameters),
59 #ifdef LIBMESH_HAVE_VTK
60  _data(NULL),
61  _algorithm(NULL),
62 #endif
63  _is_pars(parameters),
64  _is_console((parameters.getCheckedPointerParam<MooseApp *>("_moose_app"))->getOutputWarehouse())
65 
66 {
67 #ifndef LIBMESH_HAVE_VTK
68  // This should be impossible to reach, the registration of MeshBaseImageSampler is also guarded
69  // with LIBMESH_HAVE_VTK
70  mooseError("libMesh must be configured with VTK enabled to utilize MeshBaseImageSampler");
71 #endif
72 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
T getCheckedPointerParam(const std::string &name, const std::string &error_string="") const
Verifies that the requested parameter exists and is not NULL and returns it to the caller...
Base class for MOOSE-based applications.
Definition: MooseApp.h:58
vtkImageData * _data
Complete image data.
FileRangeBuilder(const InputParameters &params)
ConsoleStream _is_console
Create a console stream object for this helper class.

Member Function Documentation

◆ errorCheck()

void FileRangeBuilder::errorCheck ( )
protectedinherited

Definition at line 159 of file FileRangeBuilder.C.

Referenced by ImageMesh::ImageMesh().

160 {
161  switch (_status)
162  {
163  case 0:
164  return;
165  case 1:
166  mooseError("Cannot provide both file and file_base parameters");
167  break;
168  case 2:
169  mooseError("You must provide a valid value for either the 'file' parameter or the "
170  "'file_base' parameter.");
171  break;
172  case 3:
173  mooseError(
174  "If you provide a 'file_base', you must also provide a valid 'file_suffix', e.g. 'png'.");
175  break;
176  default:
177  mooseError("Unknown error code!");
178  }
179 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208

◆ filenames()

const std::vector<std::string>& FileRangeBuilder::filenames ( )
inlineinherited

◆ fileSuffix()

std::string FileRangeBuilder::fileSuffix ( )
inlineinherited

Definition at line 45 of file FileRangeBuilder.h.

Referenced by ImageSampler::setupImageSampler(), and setupImageSampler().

45 { return _file_suffix; }
std::string _file_suffix

◆ imageFlip()

vtkSmartPointer< vtkImageFlip > MeshBaseImageSampler::imageFlip ( const int &  axis)
private

Helper method for flipping image.

Parameters
axisFlag for determing the flip axis: "x=0", "y=1", "z=2"
Returns
A smart pointer the flipping filter

Definition at line 357 of file MeshBaseImageSampler.C.

Referenced by vtkFlip().

358 {
359  vtkSmartPointer<vtkImageFlip> flip_image = vtkSmartPointer<vtkImageFlip>::New();
360 
361  flip_image->SetFilteredAxis(axis);
362 
363  // Set the data source
364  flip_image->SetInputConnection(_algorithm);
365 
366  // Perform the flip
367  flip_image->Update();
368 
369  // Return the flip filter pointer
370  return flip_image;
371 }
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...

◆ sample()

Real MeshBaseImageSampler::sample ( const Point &  p)
virtual

Return the pixel value for the given point.

Parameters
pThe point at which to extract pixel data

Definition at line 216 of file MeshBaseImageSampler.C.

Referenced by ImageSubdomainGenerator::generate().

217 {
218 #ifdef LIBMESH_HAVE_VTK
219 
220  // Do nothing if the point is outside of the image domain
221  if (!_bounding_box.contains_point(p))
222  return 0.0;
223 
224  // Determine pixel coordinates
225  std::vector<int> x(3, 0);
226  for (int i = 0; i < LIBMESH_DIM; ++i)
227  {
228  // Compute position, only if voxel size is greater than zero
229  if (_voxel[i] == 0)
230  x[i] = 0;
231 
232  else
233  {
234  x[i] = std::floor((p(i) - _origin(i)) / _voxel[i]);
235 
236  // If the point falls on the mesh extents the index needs to be decreased by one
237  if (x[i] == _dims[i])
238  x[i]--;
239  }
240  }
241 
242  // Return the image data at the given point
243  return _data->GetScalarComponentAsDouble(x[0], x[1], x[2], _component);
244 
245 #else
246  libmesh_ignore(p); // avoid un-used parameter warnings
247  return 0.0;
248 #endif
249 }
std::vector< int > _dims
Pixel dimension of image.
Point _origin
Origin of image.
static PetscErrorCode Vec x
BoundingBox _bounding_box
Bounding box for testing points.
unsigned int _component
Component to extract.
vtkImageData * _data
Complete image data.
std::vector< double > _voxel
Physical pixel size.

◆ setupImageSampler()

void MeshBaseImageSampler::setupImageSampler ( MeshBase &  mesh)
virtual

Perform initialization of image data.

Definition at line 75 of file MeshBaseImageSampler.C.

Referenced by ImageSubdomainGenerator::generate().

76 {
77  // Don't warn that mesh or _is_pars are unused when VTK is not enabled.
78  libmesh_ignore(mesh);
79  libmesh_ignore(_is_pars);
80 
81 #ifdef LIBMESH_HAVE_VTK
82  // Get access to the Mesh object
83  BoundingBox bbox = MeshTools::create_bounding_box(mesh);
84 
85  // Set the dimensions from the Mesh if not set by the User
86  if (_is_pars.isParamValid("dimensions"))
87  _physical_dims = _is_pars.get<Point>("dimensions");
88 
89  else
90  {
91  _physical_dims(0) = bbox.max()(0) - bbox.min()(0);
92 #if LIBMESH_DIM > 1
93  _physical_dims(1) = bbox.max()(1) - bbox.min()(1);
94 #endif
95 #if LIBMESH_DIM > 2
96  _physical_dims(2) = bbox.max()(2) - bbox.min()(2);
97 #endif
98  }
99 
100  // Set the origin from the Mesh if not set in the input file
101  if (_is_pars.isParamValid("origin"))
102  _origin = _is_pars.get<Point>("origin");
103  else
104  {
105  _origin(0) = bbox.min()(0);
106 #if LIBMESH_DIM > 1
107  _origin(1) = bbox.min()(1);
108 #endif
109 #if LIBMESH_DIM > 2
110  _origin(2) = bbox.min()(2);
111 #endif
112  }
113 
114  // An array of filenames, to be filled in
115  std::vector<std::string> filenames;
116 
117  // The file suffix, to be determined
118  std::string file_suffix;
119 
120  // Try to parse our own file range parameters. If that fails, then
121  // see if the associated Mesh is an ImageMesh and use its. If that
122  // also fails, then we have to throw an error...
123  //
124  // The parseFileRange method sets parameters, thus a writable reference to the InputParameters
125  // object must be obtained from the warehouse. Generally, this should be avoided, but
126  // this is a special case.
127  if (_status != 0)
128  {
129  // TO DO : enable this. It was taken from ImageSampler.C, but cn't be implemented
130  // the same way.
131  // We don't have parameters, so see if we can get them from ImageMesh
132  /*ImageMeshgenerator * image_mesh_gen = dynamic_cast<ImageMesh *>(&mesh);
133  if (!image_mesh)
134  mooseError("No file range parameters were provided and the Mesh is not an ImageMesh.");
135 
136  // Get the ImageMesh's parameters. This should work, otherwise
137  // errors would already have been thrown...
138  filenames = image_mesh->filenames();
139  file_suffix = image_mesh->fileSuffix();*/
140  }
141  else
142  {
143  // Use our own parameters (using 'this' b/c of conflicts with filenames the local variable)
144  filenames = this->filenames();
145  file_suffix = fileSuffix();
146  }
147 
148  // Storage for the file names
149  _files = vtkSmartPointer<vtkStringArray>::New();
150 
151  for (const auto & filename : filenames)
152  _files->InsertNextValue(filename);
153 
154  // Error if no files where located
155  if (_files->GetNumberOfValues() == 0)
156  mooseError("No image file(s) located");
157 
158  // Read the image stack. Hurray for VTK not using polymorphism in a
159  // smart way... we actually have to explicitly create the type of
160  // reader based on the file extension, using an if-statement...
161  if (file_suffix == "png")
162  _image = vtkSmartPointer<vtkPNGReader>::New();
163  else if (file_suffix == "tiff" || file_suffix == "tif")
164  _image = vtkSmartPointer<vtkTIFFReader>::New();
165  else
166  mooseError("Un-supported file type '", file_suffix, "'");
167 
168  // Now that _image is set up, actually read the images
169  // Indicate that data read has started
170  _is_console << "Reading image(s)..." << std::endl;
171 
172  // Extract the data
173  _image->SetFileNames(_files);
174  _image->Update();
175  _data = _image->GetOutput();
176  _algorithm = _image->GetOutputPort();
177 
178  // Set the image dimensions and voxel size member variable
179  int * dims = _data->GetDimensions();
180  for (unsigned int i = 0; i < 3; ++i)
181  {
182  _dims.push_back(dims[i]);
183  _voxel.push_back(_physical_dims(i) / _dims[i]);
184  }
185 
186  // Set the dimensions of the image and bounding box
187  _data->SetSpacing(_voxel[0], _voxel[1], _voxel[2]);
188  _data->SetOrigin(_origin(0), _origin(1), _origin(2));
189  _bounding_box.min() = _origin;
191 
192  // Indicate data read is completed
193  _is_console << " ...image read finished" << std::endl;
194 
195  // Set the component parameter
196  // If the parameter is not set then vtkMagnitude() will applied
197  if (_is_pars.isParamValid("component"))
198  {
199  unsigned int n = _data->GetNumberOfScalarComponents();
200  _component = _is_pars.get<unsigned int>("component");
201  if (_component >= n)
202  mooseError("'component' parameter must be empty or have a value of 0 to ", n - 1);
203  }
204  else
205  _component = 0;
206 
207  // Apply filters, the toggling on and off of each filter is handled internally
208  vtkMagnitude();
210  vtkThreshold();
211  vtkFlip();
212 #endif
213 }
void vtkMagnitude()
Convert the image to greyscale.
Point _physical_dims
Physical dimensions of image.
std::vector< int > _dims
Pixel dimension of image.
Point _origin
Origin of image.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
const std::vector< std::string > & filenames()
BoundingBox _bounding_box
Bounding box for testing points.
std::string fileSuffix()
unsigned int _component
Component to extract.
vtkSmartPointer< vtkStringArray > _files
List of file names to extract data.
PetscInt n
vtkImageData * _data
Complete image data.
std::vector< double > _voxel
Physical pixel size.
void vtkFlip()
Perform image flipping.
void vtkThreshold()
Perform thresholding.
vtkSmartPointer< vtkImageReader2 > _image
Complete image data.
ConsoleStream _is_console
Create a console stream object for this helper class.
void vtkShiftAndScale()
Apply image re-scaling using the vtkImageShiftAndRescale object.
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.

◆ vtkFlip()

void MeshBaseImageSampler::vtkFlip ( )
protected

Perform image flipping.

Flip the image along the x, y, and/or z axis. If multiple flips occur, they happen in order.

Definition at line 334 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

335 {
336 #ifdef LIBMESH_HAVE_VTK
337  // Convert boolean values into an integer array, then loop over it
338  int mask[3] = {
339  _is_pars.get<bool>("flip_x"), _is_pars.get<bool>("flip_y"), _is_pars.get<bool>("flip_z")};
340 
341  for (int dim = 0; dim < 3; ++dim)
342  {
343  if (mask[dim])
344  {
345  _flip_filter = imageFlip(dim);
346 
347  // Update pointers
348  _data = _flip_filter->GetOutput();
349  _algorithm = _flip_filter->GetOutputPort();
350  }
351  }
352 #endif
353 }
vtkSmartPointer< vtkImageFlip > imageFlip(const int &axis)
Helper method for flipping image.
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
vtkSmartPointer< vtkImageFlip > _flip_filter
Pointers to image flipping filter. May be used for x, y, or z.
vtkImageData * _data
Complete image data.

◆ vtkMagnitude()

void MeshBaseImageSampler::vtkMagnitude ( )
protected

Convert the image to greyscale.

By leaving the 'component' input parameter empty, this is called automatically.

Definition at line 252 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

253 {
254 #ifdef LIBMESH_HAVE_VTK
255  // Do nothing if 'component' is set
256  if (_is_pars.isParamValid("component"))
257  return;
258 
259  // Apply the greyscale filtering
260  _magnitude_filter = vtkSmartPointer<vtkImageMagnitude>::New();
261  _magnitude_filter->SetInputConnection(_algorithm);
262  _magnitude_filter->Update();
263 
264  // Update the pointers
265  _data = _magnitude_filter->GetOutput();
266  _algorithm = _magnitude_filter->GetOutputPort();
267 #endif
268 }
vtkSmartPointer< vtkImageMagnitude > _magnitude_filter
Pointer to the magnitude filter.
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
vtkImageData * _data
Complete image data.
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.

◆ vtkShiftAndScale()

void MeshBaseImageSampler::vtkShiftAndScale ( )
protected

Apply image re-scaling using the vtkImageShiftAndRescale object.

Definition at line 271 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

272 {
273 #ifdef LIBMESH_HAVE_VTK
274  // Capture the parameters
275  double shift = _is_pars.get<double>("shift");
276  double scale = _is_pars.get<double>("scale");
277 
278  // Do nothing if shift and scale are not set
279  if (shift == 0 && scale == 1)
280  return;
281 
282  // Perform the scaling and offset actions
283  _shift_scale_filter = vtkSmartPointer<vtkImageShiftScale>::New();
284  _shift_scale_filter->SetOutputScalarTypeToDouble();
285 
286  _shift_scale_filter->SetInputConnection(_algorithm);
287  _shift_scale_filter->SetShift(shift);
288  _shift_scale_filter->SetScale(scale);
289  _shift_scale_filter->Update();
290 
291  // Update the pointers
292  _data = _shift_scale_filter->GetOutput();
293  _algorithm = _shift_scale_filter->GetOutputPort();
294 #endif
295 }
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
vtkImageData * _data
Complete image data.
vtkSmartPointer< vtkImageShiftScale > _shift_scale_filter
Pointer to the shift and scaling filter.

◆ vtkThreshold()

void MeshBaseImageSampler::vtkThreshold ( )
protected

Perform thresholding.

Definition at line 298 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

299 {
300 #ifdef LIBMESH_HAVE_VTK
301  // Do nothing if threshold not set
302  if (!_is_pars.isParamValid("threshold"))
303  return;
304 
305  // Error if both upper and lower are not set
306  if (!_is_pars.isParamValid("upper_value") || !_is_pars.isParamValid("lower_value"))
307  mooseError("When thresholding is applied, both the upper_value and lower_value parameters must "
308  "be set");
309 
310  // Create the thresholding object
311  _image_threshold = vtkSmartPointer<vtkImageThreshold>::New();
312 
313  // Set the data source
314  _image_threshold->SetInputConnection(_algorithm);
315 
316  // Setup the thresholding options
317  _image_threshold->ThresholdByUpper(_is_pars.get<Real>("threshold"));
318  _image_threshold->ReplaceInOn();
319  _image_threshold->SetInValue(_is_pars.get<Real>("upper_value"));
320  _image_threshold->ReplaceOutOn();
321  _image_threshold->SetOutValue(_is_pars.get<Real>("lower_value"));
322  _image_threshold->SetOutputScalarTypeToDouble();
323 
324  // Perform the thresholding
325  _image_threshold->Update();
326 
327  // Update the pointers
328  _data = _image_threshold->GetOutput();
329  _algorithm = _image_threshold->GetOutputPort();
330 #endif
331 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
vtkImageData * _data
Complete image data.
vtkSmartPointer< vtkImageThreshold > _image_threshold
Pointer to thresholding filter.
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.

Member Data Documentation

◆ _algorithm

vtkAlgorithmOutput* MeshBaseImageSampler::_algorithm
private

VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...

Definition at line 111 of file MeshBaseImageSampler.h.

Referenced by imageFlip(), setupImageSampler(), vtkFlip(), vtkMagnitude(), vtkShiftAndScale(), and vtkThreshold().

◆ _bounding_box

BoundingBox MeshBaseImageSampler::_bounding_box
private

Bounding box for testing points.

Definition at line 156 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _component

unsigned int MeshBaseImageSampler::_component
private

Component to extract.

Definition at line 152 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _data

vtkImageData* MeshBaseImageSampler::_data
private

Complete image data.

Definition at line 108 of file MeshBaseImageSampler.h.

Referenced by sample(), setupImageSampler(), vtkFlip(), vtkMagnitude(), vtkShiftAndScale(), and vtkThreshold().

◆ _dims

std::vector<int> MeshBaseImageSampler::_dims
private

Pixel dimension of image.

Definition at line 142 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _file_suffix

std::string FileRangeBuilder::_file_suffix
protectedinherited

◆ _filenames

std::vector<std::string> FileRangeBuilder::_filenames
protectedinherited

◆ _files

vtkSmartPointer<vtkStringArray> MeshBaseImageSampler::_files
private

List of file names to extract data.

Definition at line 105 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _flip_filter

vtkSmartPointer<vtkImageFlip> MeshBaseImageSampler::_flip_filter
private

Pointers to image flipping filter. May be used for x, y, or z.

Definition at line 126 of file MeshBaseImageSampler.h.

Referenced by vtkFlip().

◆ _image

vtkSmartPointer<vtkImageReader2> MeshBaseImageSampler::_image
private

Complete image data.

Definition at line 114 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _image_threshold

vtkSmartPointer<vtkImageThreshold> MeshBaseImageSampler::_image_threshold
private

Pointer to thresholding filter.

Definition at line 117 of file MeshBaseImageSampler.h.

Referenced by vtkThreshold().

◆ _is_console

ConsoleStream MeshBaseImageSampler::_is_console
private

Create a console stream object for this helper class.

Definition at line 162 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _is_pars

const InputParameters& MeshBaseImageSampler::_is_pars
private

Parameters for interface.

Definition at line 159 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler(), vtkFlip(), vtkMagnitude(), vtkShiftAndScale(), and vtkThreshold().

◆ _magnitude_filter

vtkSmartPointer<vtkImageMagnitude> MeshBaseImageSampler::_magnitude_filter
private

Pointer to the magnitude filter.

Definition at line 123 of file MeshBaseImageSampler.h.

Referenced by vtkMagnitude().

◆ _origin

Point MeshBaseImageSampler::_origin
private

Origin of image.

Definition at line 139 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _physical_dims

Point MeshBaseImageSampler::_physical_dims
private

Physical dimensions of image.

Definition at line 145 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _shift_scale_filter

vtkSmartPointer<vtkImageShiftScale> MeshBaseImageSampler::_shift_scale_filter
private

Pointer to the shift and scaling filter.

Definition at line 120 of file MeshBaseImageSampler.h.

Referenced by vtkShiftAndScale().

◆ _status

int FileRangeBuilder::_status
protectedinherited

◆ _voxel

std::vector<double> MeshBaseImageSampler::_voxel
private

Physical pixel size.

Definition at line 148 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().


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