ITK  5.4.0
Insight Toolkit
itkSpatialFunctionImageEvaluatorFilter.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18#ifndef itkSpatialFunctionImageEvaluatorFilter_h
19#define itkSpatialFunctionImageEvaluatorFilter_h
20
21#include "itkImageFunction.h"
24#include "itkSize.h"
25#include "itkSpatialFunction.h"
26
27namespace itk
28{
42template <typename TSpatialFunction, typename TInputImage, typename TOutputImage>
43class ITK_TEMPLATE_EXPORT SpatialFunctionImageEvaluatorFilter : public ImageToImageFilter<TInputImage, TOutputImage>
44{
45public:
46 ITK_DISALLOW_COPY_AND_MOVE(SpatialFunctionImageEvaluatorFilter);
47
53
55 itkNewMacro(Self);
56
58 itkOverrideGetNameOfClassMacro(SpatialFunctionImageEvaluatorFilter);
59
61 static constexpr unsigned int NDimensions = TInputImage::ImageDimension;
62
65
68
70 using PixelType = typename TOutputImage::PixelType;
71
74
76 using FunctionType = TSpatialFunction;
77
79 using FunctionValueType = typename FunctionType::OutputType;
80
82 using PositionType = typename FunctionType::InputType;
83
85 void
86 SetFunction(FunctionType * PixelFunction)
87 {
88 m_PixelFunction = PixelFunction;
89 }
90
91protected:
94
96 void
97 GenerateData() override;
98
99 void
100 PrintSelf(std::ostream & os, Indent indent) const override;
101
102private:
104 FunctionType * m_PixelFunction{};
105};
106} // end namespace itk
107
108#ifndef ITK_MANUAL_INSTANTIATION
109# include "itkSpatialFunctionImageEvaluatorFilter.hxx"
110#endif
111
112#endif
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Evaluates a SpatialFunction onto a source image.
~SpatialFunctionImageEvaluatorFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....