www.mooseframework.org
MeshBaseImageSampler.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 // MOOSE includes
11 #include "MeshBaseImageSampler.h"
12 #include "MooseApp.h"
13 #include "ImageMeshGenerator.h"
14 
15 #include "libmesh/mesh_tools.h"
16 
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 }
55 
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 }
72 
73 void
75 {
76  // Don't warn that mesh or _is_pars are unused when VTK is not enabled.
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));
188  _bounding_box.min() = _origin;
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 }
213 
214 Real
216 {
217 #ifdef LIBMESH_HAVE_VTK
218 
219  // Do nothing if the point is outside of the image domain
220  if (!_bounding_box.contains_point(p))
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 }
249 
250 void
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 }
268 
269 void
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 }
295 
296 void
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 }
331 
332 void
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 }
353 
354 #ifdef LIBMESH_HAVE_VTK
355 vtkSmartPointer<vtkImageFlip>
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 }
371 #endif
vtkSmartPointer< vtkImageMagnitude > _magnitude_filter
Pointer to the magnitude filter.
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.
vtkSmartPointer< vtkImageFlip > imageFlip(const int &axis)
Helper method for flipping image.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
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.
const std::vector< std::string > & filenames()
MeshBase & mesh
Base class for MOOSE-based applications.
Definition: MooseApp.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:148
static InputParameters validParams()
InputParameters emptyInputParameters()
BoundingBox _bounding_box
Bounding box for testing points.
static InputParameters validParams()
Constructor.
std::string fileSuffix()
void libmesh_ignore(const Args &...)
virtual void setupImageSampler(MeshBase &mesh)
Perform initialization of image data.
MeshBaseImageSampler(const InputParameters &parameters)
vtkSmartPointer< vtkImageFlip > _flip_filter
Pointers to image flipping filter. May be used for x, y, or z.
unsigned int _component
Component to extract.
vtkSmartPointer< vtkStringArray > _files
List of file names to extract data.
To be called in the validParams functions of classes that need to operate on ranges of files...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
vtkImageData * _data
Complete image data.
vtkSmartPointer< vtkImageShiftScale > _shift_scale_filter
Pointer to the shift and scaling filter.
vtkSmartPointer< vtkImageThreshold > _image_threshold
Pointer to thresholding filter.
std::vector< double > _voxel
Physical pixel size.
void vtkFlip()
Perform image flipping.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
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.
virtual Real sample(const Point &p)
Return the pixel value for the given point.
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...
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.