go home Home | Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itkBSplineStackTransform.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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 * http://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 itkBSplineStackTransform_h
19#define itkBSplineStackTransform_h
20
21#include "itkStackTransform.h"
23#include "elxElastixBase.h"
24
25namespace itk
26{
27template <unsigned int NDimension>
28class ITK_TEMPLATE_EXPORT BSplineStackTransform
29 : public itk::StackTransform<elx::ElastixBase::CoordRepType, NDimension, NDimension>
30{
31private:
32 using CoordRepType = elx::ElastixBase::CoordRepType;
33
34public:
36
39 using Pointer = itk::SmartPointer<BSplineStackTransform>;
40 using typename Superclass::FixedParametersType;
41 itkNewMacro(Self);
43
44private:
45 using Superclass::NumberOfGeneralFixedParametersOfStack;
46
47 static constexpr unsigned int NumberOfFixedParametersOfSubTransform =
48 AdvancedBSplineDeformableTransformBase<CoordRepType, NDimension - 1>::NumberOfFixedParameters;
49
50 static constexpr unsigned int NumberOfFixedParameters =
51 NumberOfGeneralFixedParametersOfStack + NumberOfFixedParametersOfSubTransform + 1;
52
53public:
54 void
55 SetSplineOrder(const unsigned newValue)
56 {
57 m_SplineOrder = newValue;
58
59 if (!Superclass::m_FixedParameters.empty())
60 {
61 Superclass::m_FixedParameters.back() = m_SplineOrder;
62 }
63 }
64
65protected:
68
70 ~BSplineStackTransform() override = default;
71
72 void
73 SetFixedParameters(const FixedParametersType & fixedParameters) override
74 {
75 const auto numberOfFixedParameters = fixedParameters.size();
76 if (numberOfFixedParameters < NumberOfFixedParameters)
77 {
78 itkExceptionMacro("The number of FixedParameters (" << numberOfFixedParameters << ") should be at least "
79 << NumberOfFixedParameters);
80 }
81 const auto lastFixedParameter = fixedParameters.back();
82 if (lastFixedParameter >= 1 && lastFixedParameter <= 3 &&
83 static_cast<double>(static_cast<unsigned>(lastFixedParameter)) == lastFixedParameter)
84 {
85 m_SplineOrder = static_cast<unsigned>(fixedParameters.back());
86 }
87 else
88 {
89 itkExceptionMacro("The last FixedParameters (" << lastFixedParameter << ") should be a valid spline order.");
90 }
91
92 if (Superclass::m_FixedParameters != fixedParameters)
93 {
94 Superclass::m_FixedParameters = fixedParameters;
95
96 Superclass::CreateSubTransforms(FixedParametersType(
97 fixedParameters.data_block() + NumberOfGeneralFixedParametersOfStack, NumberOfFixedParametersOfSubTransform));
98 Superclass::UpdateStackSpacingAndOrigin();
99 this->Modified();
100 }
101 }
102
103private:
104 void
105 UpdateFixedParametersInternally(const FixedParametersType & fixedParametersOfSubTransform) override
106 {
107 FixedParametersType & fixedParametersOfStack = this->Superclass::m_FixedParameters;
108 fixedParametersOfStack.set_size(NumberOfGeneralFixedParametersOfStack + NumberOfFixedParametersOfSubTransform + 1);
109 fixedParametersOfStack.back() = m_SplineOrder;
110 Superclass::UpdateFixedParametersInternally(fixedParametersOfSubTransform);
111 }
112
114 typename Superclass::SubTransformPointer
115 CreateSubTransform() const override
116 {
117 return AdvancedBSplineDeformableTransformBase<CoordRepType, NDimension - 1>::template Create<
119 }
120
121 unsigned m_SplineOrder{};
122};
123
124} // namespace itk
125
126#endif
Base class for deformable transform using a B-spline representation.
Deformable transform using a B-spline representation.
void SetSplineOrder(const unsigned newValue)
itk::SmartPointer< BSplineStackTransform > Pointer
elx::ElastixBase::CoordRepType CoordRepType
~BSplineStackTransform() override=default
void UpdateFixedParametersInternally(const FixedParametersType &fixedParametersOfSubTransform) override
void SetFixedParameters(const FixedParametersType &fixedParameters) override
Superclass::SubTransformPointer CreateSubTransform() const override
ITK_DISALLOW_COPY_AND_MOVE(BSplineStackTransform)
Implements stack of transforms: one for every last dimension index.


Generated on 2024-07-17 for elastix by doxygen 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) elastix logo