ITK  6.0.0
Insight Toolkit
itkCommandIterationUpdate.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 itkCommandIterationUpdate_h
19 #define itkCommandIterationUpdate_h
20 
21 #include "itkCommand.h"
22 #include "itkWeakPointer.h"
23 
24 namespace itk
25 {
26 
32 template <typename TOptimizer>
34 {
35 public:
36 
41 
42 
47 
48 
53 
58 
62  void
63  Execute(itk::Object * caller, const itk::EventObject & event) override
64  {
65  Execute((const itk::Object *)caller, event);
66  }
67 
68  void
69  Execute(const itk::Object *, const itk::EventObject & event) override
70  {
71  if (typeid(event) == typeid(itk::StartEvent))
72  {
73  std::cout << std::endl << "Position Value" << std::endl << std::endl;
74  }
75  else if (typeid(event) == typeid(itk::IterationEvent))
76  {
77  std::cout << m_Optimizer->GetCurrentIteration() << " = " << m_Optimizer->GetValue() << " : "
78  << m_Optimizer->GetCurrentPosition() << std::endl;
79  }
80  else if (typeid(event) == typeid(itk::EndEvent))
81  {
82  std::cout << std::endl
83  << std::endl
84  << "After " << m_Optimizer->GetCurrentIteration() << " iterations " << std::endl
85  << "Solution is = " << m_Optimizer->GetCurrentPosition() << std::endl
86  << "With value = " << m_Optimizer->GetValue() << std::endl
87  << "Stop condition = " << m_Optimizer->GetStopCondition() << std::endl;
88  }
89  }
90 
91 
95  itkOverrideGetNameOfClassMacro(CommandIterationUpdate);
96 
97 
101  itkNewMacro(Self);
102 
103 
107  using OptimizerType = TOptimizer;
108 
109 
113  void
115  {
116  m_Optimizer = optimizer;
117  m_Optimizer->AddObserver(itk::IterationEvent(), this);
118  }
121 protected:
125  CommandIterationUpdate() = default;
126 
127 private:
132 };
133 
134 } // end namespace itk
135 
136 #endif
itk::CommandIterationUpdate::m_Optimizer
WeakPointer< OptimizerType > m_Optimizer
Definition: itkCommandIterationUpdate.h:131
itkWeakPointer.h
itk::CommandIterationUpdate::Execute
void Execute(itk::Object *caller, const itk::EventObject &event) override
Definition: itkCommandIterationUpdate.h:63
itk::SmartPointer< Self >
itk::CommandIterationUpdate::CommandIterationUpdate
CommandIterationUpdate()=default
itk::Command
Superclass for callback/observer methods.
Definition: itkCommand.h:45
itk::Command
class ITK_FORWARD_EXPORT Command
Definition: itkObject.h:42
itk::CommandIterationUpdate
Definition: itkCommandIterationUpdate.h:33
itk::WeakPointer< OptimizerType >
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::EventObject
Abstraction of the Events used to communicating among filters and with GUIs.
Definition: itkEventObject.h:58
itk::CommandIterationUpdate::Execute
void Execute(const itk::Object *, const itk::EventObject &event) override
Definition: itkCommandIterationUpdate.h:69
AddImageFilter
Definition: itkAddImageFilter.h:81
itkCommand.h
itk::CommandIterationUpdate::OptimizerType
TOptimizer OptimizerType
Definition: itkCommandIterationUpdate.h:107
itk::CommandIterationUpdate::SetOptimizer
void SetOptimizer(OptimizerType *optimizer)
Definition: itkCommandIterationUpdate.h:114