RTK  2.5.0
Reconstruction Toolkit
rtkIterationCommands.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK 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  * 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 
19 #ifndef rtkIterationCommands_h
20 #define rtkIterationCommands_h
21 
22 #include <itkImageFileWriter.h>
23 
24 namespace rtk
25 {
26 
37 template <typename TCaller>
38 class ITK_TEMPLATE_EXPORT IterationCommand : public itk::Command
39 {
40 public:
45 
46  itkSetMacro(TriggerEvery, unsigned int);
47 
48  void
49  Execute(itk::Object * caller, const itk::EventObject & event) override
50  {
51  Execute((const itk::Object *)caller, event);
52  }
53 
54  void
55  Execute(const itk::Object * caller, const itk::EventObject & event) override
56  {
57  const auto * cCaller = dynamic_cast<const TCaller *>(caller);
58  if (cCaller)
59  {
60  if (itk::EndEvent().CheckEvent(&event))
61  {
62  End(cCaller);
63  return;
64  }
65  if (itk::IterationEvent().CheckEvent(&event))
66  {
67  ++m_IterationCount;
68  if ((m_IterationCount % m_TriggerEvery) == 0)
69  {
70  Run(cCaller);
71  }
72  }
73  } // TODO fail when cast fails
74  }
75 
76 protected:
78  unsigned int m_IterationCount = 0;
79 
81  unsigned int m_TriggerEvery = 1;
82 
84  virtual void
85  Run(const TCaller * caller) = 0;
86 
88  virtual void
89  End(const TCaller * itkNotUsed(caller))
90  { /* Default implementation: do nothing */
91  }
92 };
93 
102 template <typename TCaller>
103 class ITK_TEMPLATE_EXPORT VerboseIterationCommand : public IterationCommand<TCaller>
104 {
105 public:
110  itkNewMacro(Self);
111 
112 protected:
113  void
114  Run(const TCaller * itkNotUsed(caller)) override
115  {
116  // TODO allow string personnalization?
117  std::cout << "\rIteration " << this->m_IterationCount << " completed." << std::flush;
118  }
119 
120  void
121  End(const TCaller * itkNotUsed(caller)) override
122  {
123  std::cout << std::endl; // new line when execution ends
124  }
125 };
126 
138 template <typename TCaller, typename TOutputImage>
139 class ITK_TEMPLATE_EXPORT OutputIterationCommand : public IterationCommand<TCaller>
140 {
141 public:
146  itkNewMacro(Self);
147 
148  itkSetMacro(FileFormat, std::string);
149 
150 protected:
152  std::string m_FileFormat;
153 
154  void
155  Run(const TCaller * caller) override
156  {
157  typedef itk::ImageFileWriter<TOutputImage> WriterType;
158  typename WriterType::Pointer writer = WriterType::New();
159 
160  char buffer[100];
161  unsigned int size = sprintf(buffer, m_FileFormat.c_str(), this->m_IterationCount);
162  writer->SetFileName(std::string(buffer, size));
163 
164  writer->SetInput(caller->GetOutput());
165  writer->Update();
166  }
167 };
168 
169 } // end namespace rtk
170 
171 #endif
IterationCommand< TCaller > Superclass
IterationCommand< TCaller > Superclass
Outputs to standard output when an iteration completes.
void Execute(const itk::Object *caller, const itk::EventObject &event) override
virtual void End(const TCaller *)
itk::SmartPointer< Self > Pointer
void End(const TCaller *) override
VerboseIterationCommand Self
Output intermediate iterations in a file. This class is useful to check convergence of an iterative m...
void Run(const TCaller *caller) override
#define itkSetMacro(name, type)
itk::SmartPointer< Self > Pointer
Abstract superclass to all iteration callbacks. Derived classes must implement the Run() method...
void Execute(itk::Object *caller, const itk::EventObject &event) override
itk::SmartPointer< Self > Pointer
void Run(const TCaller *) override