ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMultiGradientOptimizerv4.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 itkMultiGradientOptimizerv4_h
19 #define itkMultiGradientOptimizerv4_h
20 
23 
24 namespace itk
25 {
45 template<typename TInternalComputationValueType>
46 class ITK_TEMPLATE_EXPORT MultiGradientOptimizerv4Template
47 : public GradientDescentOptimizerv4Template<TInternalComputationValueType>
48 {
49 public:
50  ITK_DISALLOW_COPY_AND_ASSIGN(MultiGradientOptimizerv4Template);
51 
57 
60 
62  itkNewMacro(Self);
63 
66  using ParametersType = typename Superclass::ParametersType;
69  using OptimizersListType = std::vector< LocalOptimizerPointer >;
70  using OptimizersListSizeType = typename OptimizersListType::size_type;
71 
73 
75  using StopConditionReturnStringType = typename Superclass::StopConditionReturnStringType;
76 
78  using StopConditionDescriptionType = typename Superclass::StopConditionDescriptionType;
79 
81  using InternalComputationValueType = TInternalComputationValueType;
82 
84  using MetricType = typename Superclass::MetricType;
85  using MetricTypePointer = typename MetricType::Pointer;
86 
88  using DerivativeType = typename MetricType::DerivativeType;
89 
92  using MetricValuesListType = std::vector< MeasureType >;
93 
95  const StopConditionType & GetStopCondition() const override
96  {
97  return this->m_StopCondition;
98  }
99 
101  void StartOptimization( bool doOnlyInitialization = false ) override;
102 
105  void StopOptimization() override;
106 
109  void ResumeOptimization() override;
110 
112  const StopConditionReturnStringType GetStopConditionDescription() const override;
113 
115  OptimizersListType & GetOptimizersList();
116 
118  void SetOptimizersList(OptimizersListType & p);
119 
121  const MetricValuesListType & GetMetricValuesList() const;
122 
123  protected:
124 
127  ~MultiGradientOptimizerv4Template() override = default;
129 
130  void PrintSelf(std::ostream & os, Indent indent) const override;
131 
132  /* Common variables for optimization control and reporting */
133  bool m_Stop{false};
140 };
141 
144 
145 } // end namespace itk
146 
147 #ifndef ITK_MANUAL_INSTANTIATION
148 #include "itkMultiGradientOptimizerv4.hxx"
149 #endif
150 
151 #endif
Multiple gradient-based optimizers are combined in order to perform a multi-objective optimization...
Light weight base class for most itk classes.
std::vector< LocalOptimizerPointer > OptimizersListType
const StopConditionType & GetStopCondition() const override
typename itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >::Pointer LocalOptimizerPointer
typename OptimizersListType::size_type OptimizersListSizeType
typename Superclass::StopConditionDescriptionType StopConditionDescriptionType
typename Superclass::StopConditionReturnStringType StopConditionReturnStringType
StopConditionDescriptionType m_StopConditionDescription
typename Superclass::StopConditionType StopConditionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename OptimizerType::Pointer OptimizerPointer
Abstract base for object-to-object optimizers.