https://mooseframework.inl.gov
Public Member Functions | Static 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)
 
virtual libMesh::Real sample (const libMesh::Point &p)
 Return the pixel value for the given point. More...
 
virtual void setupImageSampler (libMesh::MeshBase &mesh)
 Perform initialization of image data. More...
 
std::string fileSuffix ()
 
const std::vector< std::string > & filenames ()
 

Static Public Member Functions

static InputParameters validParams ()
 Constructor. More...
 

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...
 
libMesh::Point _origin
 Origin of image. More...
 
std::vector< int_dims
 Pixel dimension of image. More...
 
libMesh::Point _physical_dims
 Physical dimensions of image. More...
 
std::vector< double > _voxel
 Physical pixel size. More...
 
unsigned int _component
 Component to extract. More...
 
libMesh::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 81 of file MeshBaseImageSampler.h.

Constructor & Destructor Documentation

◆ MeshBaseImageSampler()

MeshBaseImageSampler::MeshBaseImageSampler ( const InputParameters parameters)

Definition at line 56 of file MeshBaseImageSampler.C.

57  : FileRangeBuilder(parameters),
58 #ifdef LIBMESH_HAVE_VTK
59  _data(NULL),
60  _algorithm(NULL),
61 #endif
62  _is_pars(parameters),
63  _is_console((parameters.getCheckedPointerParam<MooseApp *>("_moose_app"))->getOutputWarehouse())
64 
65 {
66 #ifndef LIBMESH_HAVE_VTK
67  // This should be impossible to reach, the registration of MeshBaseImageSampler is also guarded
68  // with LIBMESH_HAVE_VTK
69  mooseError("libMesh must be configured with VTK enabled to utilize MeshBaseImageSampler");
70 #endif
71 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
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:96
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:302

◆ filenames()

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

◆ fileSuffix()

std::string FileRangeBuilder::fileSuffix ( )
inlineinherited

Definition at line 42 of file FileRangeBuilder.h.

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

42 { 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 356 of file MeshBaseImageSampler.C.

Referenced by vtkFlip().

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

◆ sample()

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

Return the pixel value for the given point.

Parameters
pThe point at which to extract pixel data

Definition at line 215 of file MeshBaseImageSampler.C.

Referenced by ImageSubdomainGenerator::generate().

216 {
217 #ifdef LIBMESH_HAVE_VTK
218 
219  // Do nothing if the point is outside of the image domain
221  return 0.0;
222 
223  // Determine pixel coordinates
224  std::vector<int> x(3, 0);
225  for (int i = 0; i < LIBMESH_DIM; ++i)
226  {
227  // Compute position, only if voxel size is greater than zero
228  if (_voxel[i] == 0)
229  x[i] = 0;
230 
231  else
232  {
233  x[i] = std::floor((p(i) - _origin(i)) / _voxel[i]);
234 
235  // If the point falls on the mesh extents the index needs to be decreased by one
236  if (x[i] == _dims[i])
237  x[i]--;
238  }
239  }
240 
241  // Return the image data at the given point
242  return _data->GetScalarComponentAsDouble(x[0], x[1], x[2], _component);
243 
244 #else
245  libmesh_ignore(p); // avoid un-used parameter warnings
246  return 0.0;
247 #endif
248 }
std::vector< int > _dims
Pixel dimension of image.
bool contains_point(const Point &) const
void libmesh_ignore(const Args &...)
libMesh::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.
libMesh::Point _origin
Origin of image.

◆ setupImageSampler()

void MeshBaseImageSampler::setupImageSampler ( libMesh::MeshBase mesh)
virtual

Perform initialization of image data.

Definition at line 74 of file MeshBaseImageSampler.C.

Referenced by ImageSubdomainGenerator::generate().

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

◆ validParams()

InputParameters MeshBaseImageSampler::validParams ( )
static

Constructor.

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

See also
ImageFunction

Definition at line 18 of file MeshBaseImageSampler.C.

Referenced by ImageSubdomainGenerator::validParams().

19 {
20  // Define the general parameters
23 
24  params.addParam<Point>("origin", "Origin of the image (defaults to mesh origin)");
25  params.addParam<Point>("dimensions",
26  "x,y,z dimensions of the image (defaults to mesh dimensions)");
27  params.addParam<unsigned int>(
28  "component",
29  "The image RGB-component to return, leaving this blank will result in a greyscale value "
30  "for the image to be created. The component number is zero based, i.e. 0 returns the first "
31  "(RED) component of the image.");
32 
33  // Shift and Scale (application of these occurs prior to threshold)
34  params.addParam<double>("shift", 0, "Value to add to all pixels; occurs prior to scaling");
35  params.addParam<double>(
36  "scale", 1, "Multiplier to apply to all pixel values; occurs after shifting");
37  params.addParamNamesToGroup("shift scale", "Rescale");
38 
39  // Threshold parameters
40  params.addParam<double>("threshold", "The threshold value");
41  params.addParam<double>(
42  "upper_value", 1, "The value to set for data greater than the threshold value");
43  params.addParam<double>(
44  "lower_value", 0, "The value to set for data less than the threshold value");
45  params.addParamNamesToGroup("threshold upper_value lower_value", "Threshold");
46 
47  // Flip image
48  params.addParam<bool>("flip_x", false, "Flip the image along the x-axis");
49  params.addParam<bool>("flip_y", false, "Flip the image along the y-axis");
50  params.addParam<bool>("flip_z", false, "Flip the image along the z-axis");
51  params.addParamNamesToGroup("flip_x flip_y flip_z", "Flip");
52 
53  return params;
54 }
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
InputParameters emptyInputParameters()
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...

◆ 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 333 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

334 {
335 #ifdef LIBMESH_HAVE_VTK
336  // Convert boolean values into an integer array, then loop over it
337  int mask[3] = {
338  _is_pars.get<bool>("flip_x"), _is_pars.get<bool>("flip_y"), _is_pars.get<bool>("flip_z")};
339 
340  for (int dim = 0; dim < 3; ++dim)
341  {
342  if (mask[dim])
343  {
345 
346  // Update pointers
347  _data = _flip_filter->GetOutput();
348  _algorithm = _flip_filter->GetOutputPort();
349  }
350  }
351 #endif
352 }
vtkSmartPointer< vtkImageFlip > imageFlip(const int &axis)
Helper method for flipping image.
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
vtkAlgorithmOutput * _algorithm
VTK-6 seems to work better in terms of "algorithm outputs" rather than vtkImageData pointers...
const InputParameters & _is_pars
Parameters for interface.
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:154
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 251 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

252 {
253 #ifdef LIBMESH_HAVE_VTK
254  // Do nothing if 'component' is set
255  if (_is_pars.isParamValid("component"))
256  return;
257 
258  // Apply the greyscale filtering
259  _magnitude_filter = vtkSmartPointer<vtkImageMagnitude>::New();
260  _magnitude_filter->SetInputConnection(_algorithm);
261  _magnitude_filter->Update();
262 
263  // Update the pointers
264  _data = _magnitude_filter->GetOutput();
265  _algorithm = _magnitude_filter->GetOutputPort();
266 #endif
267 }
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 270 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

271 {
272 #ifdef LIBMESH_HAVE_VTK
273  // Capture the parameters
274  double shift = _is_pars.get<double>("shift");
275  double scale = _is_pars.get<double>("scale");
276 
277  // Do nothing if shift and scale are not set
278  if (shift == 0 && scale == 1)
279  return;
280 
281  // Perform the scaling and offset actions
282  _shift_scale_filter = vtkSmartPointer<vtkImageShiftScale>::New();
283  _shift_scale_filter->SetOutputScalarTypeToDouble();
284 
285  _shift_scale_filter->SetInputConnection(_algorithm);
286  _shift_scale_filter->SetShift(shift);
287  _shift_scale_filter->SetScale(scale);
288  _shift_scale_filter->Update();
289 
290  // Update the pointers
291  _data = _shift_scale_filter->GetOutput();
292  _algorithm = _shift_scale_filter->GetOutputPort();
293 #endif
294 }
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
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 297 of file MeshBaseImageSampler.C.

Referenced by setupImageSampler().

298 {
299 #ifdef LIBMESH_HAVE_VTK
300  // Do nothing if threshold not set
301  if (!_is_pars.isParamValid("threshold"))
302  return;
303 
304  // Error if both upper and lower are not set
305  if (!_is_pars.isParamValid("upper_value") || !_is_pars.isParamValid("lower_value"))
306  mooseError("When thresholding is applied, both the upper_value and lower_value parameters must "
307  "be set");
308 
309  // Create the thresholding object
310  _image_threshold = vtkSmartPointer<vtkImageThreshold>::New();
311 
312  // Set the data source
313  _image_threshold->SetInputConnection(_algorithm);
314 
315  // Setup the thresholding options
316  _image_threshold->ThresholdByUpper(_is_pars.get<Real>("threshold"));
317  _image_threshold->ReplaceInOn();
318  _image_threshold->SetInValue(_is_pars.get<Real>("upper_value"));
319  _image_threshold->ReplaceOutOn();
320  _image_threshold->SetOutValue(_is_pars.get<Real>("lower_value"));
321  _image_threshold->SetOutputScalarTypeToDouble();
322 
323  // Perform the thresholding
324  _image_threshold->Update();
325 
326  // Update the pointers
327  _data = _image_threshold->GetOutput();
328  _algorithm = _image_threshold->GetOutputPort();
329 #endif
330 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
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 142 of file MeshBaseImageSampler.h.

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

◆ _bounding_box

libMesh::BoundingBox MeshBaseImageSampler::_bounding_box
private

Bounding box for testing points.

Definition at line 187 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _component

unsigned int MeshBaseImageSampler::_component
private

Component to extract.

Definition at line 183 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _data

vtkImageData* MeshBaseImageSampler::_data
private

Complete image data.

Definition at line 139 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 173 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 136 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 157 of file MeshBaseImageSampler.h.

Referenced by vtkFlip().

◆ _image

vtkSmartPointer<vtkImageReader2> MeshBaseImageSampler::_image
private

Complete image data.

Definition at line 145 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _image_threshold

vtkSmartPointer<vtkImageThreshold> MeshBaseImageSampler::_image_threshold
private

Pointer to thresholding filter.

Definition at line 148 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 193 of file MeshBaseImageSampler.h.

Referenced by setupImageSampler().

◆ _is_pars

const InputParameters& MeshBaseImageSampler::_is_pars
private

Parameters for interface.

Definition at line 190 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 154 of file MeshBaseImageSampler.h.

Referenced by vtkMagnitude().

◆ _origin

libMesh::Point MeshBaseImageSampler::_origin
private

Origin of image.

Definition at line 170 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().

◆ _physical_dims

libMesh::Point MeshBaseImageSampler::_physical_dims
private

Physical dimensions of image.

Definition at line 176 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 151 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 179 of file MeshBaseImageSampler.h.

Referenced by sample(), and setupImageSampler().


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