ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkBSplineSmoothingOnUpdateDisplacementFieldTransform.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 itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
19#define itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
20
22
24#include "itkPointSet.h"
25
26namespace itk
27{
28
46template <typename TParametersValueType, unsigned int VDimension>
48 : public DisplacementFieldTransform<TParametersValueType, VDimension>
49{
50public:
52
58
61
63 itkNewMacro(Self);
64
66 static constexpr unsigned int Dimension = VDimension;
67
69 using typename Superclass::ScalarType;
70 using typename Superclass::DerivativeType;
75
77
84 using SplineOrderType = unsigned int;
89 using ArrayValueType = typename ArrayType::ValueType;
90
101 void
102 UpdateTransformParameters(const DerivativeType & update, ScalarType factor = 1.0) override;
103
107 itkSetMacro(SplineOrder, SplineOrderType);
108
112 itkGetConstMacro(SplineOrder, SplineOrderType);
113
121 itkSetMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
122
130 itkGetConstMacro(NumberOfControlPointsForTheUpdateField, ArrayType);
131
138 void
140
148 itkSetMacro(NumberOfControlPointsForTheTotalField, ArrayType);
149
157 itkGetConstMacro(NumberOfControlPointsForTheTotalField, ArrayType);
158
165 void
167
172 itkBooleanMacro(EnforceStationaryBoundary);
173 itkSetMacro(EnforceStationaryBoundary, bool);
174 itkGetConstMacro(EnforceStationaryBoundary, bool);
176protected:
179
180 void
181 PrintSelf(std::ostream & os, Indent indent) const override;
182
184 typename LightObject::Pointer
185 InternalClone() const override;
186
193
194private:
199};
200
201} // end namespace itk
202
203#ifndef ITK_MANUAL_INSTANTIATION
204# include "itkBSplineSmoothingOnUpdateDisplacementFieldTransform.hxx"
205#endif
206
207#endif // itkBSplineSmoothingOnUpdateDisplacementFieldTransform_h
ParametersValueType ValueType
Definition itkArray.h:51
DisplacementFieldPointer BSplineSmoothDisplacementField(const DisplacementFieldType *, const ArrayType &)
DisplacementFieldToBSplineImageFilter< DisplacementFieldType > BSplineFilterType
LightObject::Pointer InternalClone() const override
typename Transform< TParametersValueType, VDimension, VDimension >::Pointer TransformPointer
void PrintSelf(std::ostream &os, Indent indent) const override
void UpdateTransformParameters(const DerivativeType &update, ScalarType factor=1.0) override
DisplacementFieldTransform< TParametersValueType, VDimension > Superclass
Class which takes a dense displacement field image and/or a set of points with associated displacemen...
Array< ParametersValueType > DerivativeType
typename DisplacementFieldType::Pointer DisplacementFieldPointer
ParametersValueType ScalarType
typename DisplacementFieldType::ConstPointer DisplacementFieldConstPointer
Image< OutputVectorType, Dimension > DisplacementFieldType
Control indentation during Print() invocation.
Definition itkIndent.h:50
SmartPointer< Self > Pointer
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition itkPointSet.h:82
Implements transparent reference counting.
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....