ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkLBFGSBOptimizerv4.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 itkLBFGSBOptimizerv4_h
19#define itkLBFGSBOptimizerv4_h
20
22#include "vnl/algo/vnl_lbfgsb.h"
23#include "ITKOptimizersv4Export.h"
24
25namespace itk
26{
27/* Necessary forward declaration */
37// Forward reference because of private implementation
38class ITK_FORWARD_EXPORT LBFGSBOptimizerHelperv4;
39
56class ITKOptimizersv4_EXPORT LBFGSBOptimizerv4 : public LBFGSOptimizerBasev4<vnl_lbfgsb>
57{
58public:
59 ITK_DISALLOW_COPY_AND_MOVE(LBFGSBOptimizerv4);
60
66
70
72 itkNewMacro(Self);
73
75 itkOverrideGetNameOfClassMacro(LBFGSBOptimizerv4);
76
84
89
94
96 void
98
102 {
103 return m_InitialPosition;
104 }
105
107 void
108 StartOptimization(bool doOnlyInitialization = false) override;
109
111 void
112 SetMetric(MetricType * metric) override;
113
115 void
117
118 itkGetConstReferenceMacro(LowerBound, BoundValueType);
119
121 void
123
124 itkGetConstReferenceMacro(UpperBound, BoundValueType);
125
132 void
134
135 itkGetConstReferenceMacro(BoundSelection, BoundSelectionType);
136
143 virtual void
145
146 itkGetConstMacro(CostFunctionConvergenceFactor, double);
147
149 virtual void
151
152 itkGetConstMacro(MaximumNumberOfCorrections, unsigned int);
153
155 void
156 SetScales(const ScalesType &) override;
157
160 itkGetConstReferenceMacro(InfinityNormOfProjectedGradient, double);
161
163 bool
164 CanUseScales() const override
165 {
166 return false;
167 }
168
169protected:
172 void
173 PrintSelf(std::ostream & os, Indent indent) const override;
174
176
179
181
182private:
184
189};
190} // end namespace itk
191#endif
Array class with size defined at construction time.
Definition itkArray.h:48
Control indentation during Print() invocation.
Definition itkIndent.h:50
void SetMetric(MetricType *metric) override
void SetBoundSelection(const BoundSelectionType &value)
void SetUpperBound(const BoundValueType &value)
LBFGSOptimizerBasev4< vnl_lbfgsb > Superclass
ParametersType & GetInitialPosition()
LBFGSBOptimizerHelperv4 InternalOptimizerType
virtual void SetCostFunctionConvergenceFactor(double)
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
void PrintSelf(std::ostream &os, Indent indent) const override
BoundSelectionType m_BoundSelection
unsigned int m_MaximumNumberOfCorrections
virtual void SetMaximumNumberOfCorrections(unsigned int)
~LBFGSBOptimizerv4() override
void SetLowerBound(const BoundValueType &value)
SmartPointer< const Self > ConstPointer
void SetInitialPosition(const ParametersType &param)
void StartOptimization(bool doOnlyInitialization=false) override
SmartPointer< Self > Pointer
Superclass::ScalesType ScalesType
Superclass::MetricType MetricType
void SetScales(const ScalesType &) override
Superclass::ParametersType ParametersType
bool CanUseScales() const override
Superclass::CostFunctionAdaptorType CostFunctionAdaptorType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT LBFGSBOptimizerHelperv4