ITK  5.4.0
Insight Toolkit
itkMultiGradientOptimizerv4.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 itkMultiGradientOptimizerv4_h
19#define itkMultiGradientOptimizerv4_h
20
23
24namespace itk
25{
46template <typename TInternalComputationValueType>
47class ITK_TEMPLATE_EXPORT MultiGradientOptimizerv4Template
48 : public GradientDescentOptimizerv4Template<TInternalComputationValueType>
49{
50public:
51 ITK_DISALLOW_COPY_AND_MOVE(MultiGradientOptimizerv4Template);
52
58
60 itkOverrideGetNameOfClassMacro(MultiGradientOptimizerv4Template);
61
63 itkNewMacro(Self);
64
68 using typename Superclass::ParametersType;
71 using OptimizersListType = std::vector<LocalOptimizerPointer>;
72 using OptimizersListSizeType = typename OptimizersListType::size_type;
73
75 using typename Superclass::StopConditionReturnStringType;
76
78 using typename Superclass::StopConditionDescriptionType;
79
81 using InternalComputationValueType = TInternalComputationValueType;
82
84 using typename Superclass::MetricType;
86
89
91 using typename Superclass::MeasureType;
92 using MetricValuesListType = std::vector<MeasureType>;
93
96 GetStopCondition() const override
97 {
98 return this->m_StopCondition;
99 }
100
102 void
103 StartOptimization(bool doOnlyInitialization = false) override;
104
107 void
109
112 void
114
118
122
124 void
126
130
131protected:
137 void
138 PrintSelf(std::ostream & os, Indent indent) const override;
139
140 /* Common variables for optimization control and reporting */
141 bool m_Stop{ false };
143 StopConditionDescriptionType m_StopConditionDescription{};
144 OptimizersListType m_OptimizersList{};
145 MetricValuesListType m_MetricValuesList{};
146 MeasureType m_MinimumMetricValue{};
147 MeasureType m_MaximumMetricValue{};
148};
149
152
153} // end namespace itk
154
155#ifndef ITK_MANUAL_INSTANTIATION
156# include "itkMultiGradientOptimizerv4.hxx"
157#endif
158
159#endif
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Multiple gradient-based optimizers are combined in order to perform a multi-objective optimization.
const StopConditionObjectToObjectOptimizerEnum & GetStopCondition() const override
void PrintSelf(std::ostream &os, Indent indent) const override
typename itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >::Pointer LocalOptimizerPointer
std::vector< LocalOptimizerPointer > OptimizersListType
OptimizersListType & GetOptimizersList()
void SetOptimizersList(OptimizersListType &p)
typename OptimizersListType::size_type OptimizersListSizeType
void StartOptimization(bool doOnlyInitialization=false) override
~MultiGradientOptimizerv4Template() override=default
const StopConditionReturnStringType GetStopConditionDescription() const override
typename MetricType::DerivativeType DerivativeType
const MetricValuesListType & GetMetricValuesList() const
typename OptimizerType::Pointer OptimizerPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....