ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkTimeVaryingVelocityFieldIntegrationImageFilter.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 itkTimeVaryingVelocityFieldIntegrationImageFilter_h
19#define itkTimeVaryingVelocityFieldIntegrationImageFilter_h
20
22
24
25namespace itk
26{
52template <typename TTimeVaryingVelocityField,
53 typename TDisplacementField =
54 Image<typename TTimeVaryingVelocityField::PixelType, TTimeVaryingVelocityField::ImageDimension - 1>>
56 : public ImageToImageFilter<TTimeVaryingVelocityField, TDisplacementField>
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(TimeVaryingVelocityFieldIntegrationImageFilter);
60
65
67 itkNewMacro(Self);
68
70 itkOverrideGetNameOfClassMacro(TimeVaryingVelocityFieldIntegrationImageFilter);
71
75 static constexpr unsigned int InputImageDimension = TTimeVaryingVelocityField::ImageDimension;
76
77 static constexpr unsigned int OutputImageDimension = TDisplacementField::ImageDimension;
78
79 using TimeVaryingVelocityFieldType = TTimeVaryingVelocityField;
80 using DisplacementFieldType = TDisplacementField;
81 using DisplacementFieldPointer = typename DisplacementFieldType::Pointer;
82 using VectorType = typename DisplacementFieldType::PixelType;
83 using RealType = typename VectorType::RealValueType;
84 using ScalarType = typename VectorType::ValueType;
85 using PointType = typename DisplacementFieldType::PointType;
86 using OutputRegionType = typename DisplacementFieldType::RegionType;
87
90
93
95 itkSetObjectMacro(VelocityFieldInterpolator, VelocityFieldInterpolatorType);
96 itkGetModifiableObjectMacro(VelocityFieldInterpolator, VelocityFieldInterpolatorType);
98
103 itkSetObjectMacro(DisplacementFieldInterpolator, DisplacementFieldInterpolatorType);
104 itkGetModifiableObjectMacro(DisplacementFieldInterpolator, DisplacementFieldInterpolatorType);
106
110 itkSetObjectMacro(InitialDiffeomorphism, DisplacementFieldType);
111 itkGetModifiableObjectMacro(InitialDiffeomorphism, DisplacementFieldType);
113
118 itkSetClampMacro(LowerTimeBound, RealType, 0, 1);
119
124 itkGetConstMacro(LowerTimeBound, RealType);
125
130 itkSetClampMacro(UpperTimeBound, RealType, 0, 1);
131
136 itkGetConstMacro(UpperTimeBound, RealType);
137
142 itkSetMacro(NumberOfIntegrationSteps, unsigned int);
143
148 itkGetConstMacro(NumberOfIntegrationSteps, unsigned int);
149
158 itkSetMacro(TimeBoundsAsRates, bool);
159 itkGetConstMacro(TimeBoundsAsRates, bool);
160 itkBooleanMacro(TimeBoundsAsRates);
162
163protected:
166
167 void
168 PrintSelf(std::ostream & os, Indent indent) const override;
169
170 void
172
173 void
175
176 void
178
179
181 IntegrateVelocityAtPoint(const PointType & initialSpatialPoint, const TimeVaryingVelocityFieldType * inputField);
182
185
187
189
190 unsigned int m_NumberOfTimePoints{};
191
193
195
196private:
198};
199} // namespace itk
200
201#ifndef ITK_MANUAL_INSTANTIATION
202# include "itkTimeVaryingVelocityFieldIntegrationImageFilter.hxx"
203#endif
204
205#endif
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
VectorInterpolateImageFunction< TimeVaryingVelocityFieldType, ScalarType > VelocityFieldInterpolatorType
ImageToImageFilter< TTimeVaryingVelocityField, TDisplacementField > Superclass
typename DisplacementFieldInterpolatorType::Pointer DisplacementFieldInterpolatorPointer
void PrintSelf(std::ostream &os, Indent indent) const override
void DynamicThreadedGenerateData(const OutputRegionType &) override
VectorInterpolateImageFunction< DisplacementFieldType, ScalarType > DisplacementFieldInterpolatorType
VectorType IntegrateVelocityAtPoint(const PointType &initialSpatialPoint, const TimeVaryingVelocityFieldType *inputField)
Base class for all vector image interpolators.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....