go home Home | Main Page | Topics | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
Loading...
Searching...
No Matches
elxImpactMetric.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright UMC Utrecht and contributors
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
19#ifndef elxImpactMetric_h
20#define elxImpactMetric_h
21
22#include "elxIncludes.h" // include first to avoid MSVS warning
23#include "elxConfiguration.h"
25#include "itkVector.h"
26
27namespace elastix
28{
29
177
178template <typename TElastix>
179class ITK_TEMPLATE_EXPORT ImpactMetric
180 : public itk::ImpactImageToImageMetric<typename MetricBase<TElastix>::FixedImageType,
181 typename MetricBase<TElastix>::MovingImageType>
182 , public MetricBase<TElastix>
183{
184public:
186
192 using Pointer = itk::SmartPointer<Self>;
193 using ConstPointer = itk::SmartPointer<const Self>;
194
196 itkNewMacro(Self);
197
200
210
212 using typename Superclass1::CoordinateRepresentationType;
213 using typename Superclass1::MovingImageType;
214 using typename Superclass1::MovingImagePixelType;
215 using typename Superclass1::MovingImageConstPointer;
216 using typename Superclass1::FixedImageType;
217 using typename Superclass1::FixedImageConstPointer;
218 using typename Superclass1::FixedImageRegionType;
219 using typename Superclass1::TransformType;
220 using typename Superclass1::TransformPointer;
221 using typename Superclass1::InputPointType;
222 using typename Superclass1::OutputPointType;
223 using typename Superclass1::TransformJacobianType;
224 using typename Superclass1::InterpolatorType;
225 using typename Superclass1::InterpolatorPointer;
226 using typename Superclass1::RealType;
227 using typename Superclass1::GradientPixelType;
228 using typename Superclass1::GradientImageType;
229 using typename Superclass1::GradientImagePointer;
230 using typename Superclass1::FixedImageMaskType;
231 using typename Superclass1::FixedImageMaskPointer;
232 using typename Superclass1::MovingImageMaskType;
233 using typename Superclass1::MovingImageMaskPointer;
234 using typename Superclass1::MeasureType;
235 using typename Superclass1::DerivativeType;
236 using typename Superclass1::ParametersType;
237 using typename Superclass1::FixedImagePixelType;
238 using typename Superclass1::MovingImageRegionType;
239 using typename Superclass1::ImageSamplerType;
240 using typename Superclass1::ImageSamplerPointer;
241 using typename Superclass1::ImageSampleContainerType;
242 using typename Superclass1::ImageSampleContainerPointer;
243 using typename Superclass1::FixedImageLimiterType;
244 using typename Superclass1::MovingImageLimiterType;
245 using typename Superclass1::FixedImageLimiterOutputType;
246 using typename Superclass1::MovingImageLimiterOutputType;
247 using typename Superclass1::MovingImageDerivativeScalesType;
248
250 itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension);
251
253 itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension);
254
256 using typename Superclass2::ElastixType;
257 using typename Superclass2::RegistrationType;
259
260
264 void
265 Initialize() override;
266
269 void
271
277 void
279
281 void
283
284protected:
286 ImpactMetric() = default;
288 ~ImpactMetric() override = default;
289
290 unsigned long m_CurrentIteration;
291
292private:
294
295 std::vector<itk::ImpactModelConfiguration>
296 GenerateModelsConfiguration(unsigned int level,
297 std::string type,
298 std::string mode,
299 unsigned int imageDimension,
300 bool useMixedPrecision);
301};
302
316inline std::vector<bool>
317GetBooleanVectorFromString(std::string valueStr, bool defaultValue)
318{
319 std::vector<bool> values(defaultValue, valueStr.size());
320 for (char c : valueStr)
321 {
322 std::stringstream tokenStream(std::string(1, c));
323 std::string subToken;
324 if (tokenStream >> subToken)
325 {
327 bool value = subToken == "1";
328 values.push_back(value);
329 }
330 }
331 if (values.empty())
332 {
333 values.push_back(defaultValue);
334 }
335 return values;
336} // end GetVectorFromString
337
355template <typename T>
356std::vector<std::vector<T>>
357GroupByDimensions(const std::vector<T> & values, const std::vector<unsigned int> & dimensions)
358{
359 std::vector<std::vector<T>> grouped;
360 size_t currentIndex = 0;
361 size_t n = values.size();
362
363 for (unsigned int dim : dimensions)
364 {
365 std::vector<T> group;
366 for (unsigned int i = 0; i < dim; ++i)
367 {
368 if (currentIndex < n)
369 {
370 group.push_back(values[currentIndex++]);
371 }
372 else if (!values.empty())
373 {
374 group.push_back(values.back()); // repeats last value if overflow
375 }
376 }
377 grouped.push_back(std::move(group));
378 }
379
380 return grouped;
381}
382
383} // end namespace elastix
384
385
386#ifndef ITK_MANUAL_INSTANTIATION
387# include "elxImpactMetric.hxx"
388#endif
389
390#endif // end #ifndef elxImpactMetric_h
void Initialize() override
itk::ImpactImageToImageMetric< typename MetricBase< TElastix >::FixedImageType, typename MetricBase< TElastix >::MovingImageType > Superclass1
elxClassNameMacro("Impact")
ITK_DISALLOW_COPY_AND_MOVE(ImpactMetric)
~ImpactMetric() override=default
itk::SmartPointer< Self > Pointer
itk::SmartPointer< const Self > ConstPointer
MetricBase< TElastix > Superclass2
void AfterEachIteration() override
typename Superclass2::ITKBaseType ITKBaseType
void BeforeEachResolution() override
std::vector< itk::ImpactModelConfiguration > GenerateModelsConfiguration(unsigned int level, std::string type, std::string mode, unsigned int imageDimension, bool useMixedPrecision)
void BeforeRegistration() override
itkOverrideGetNameOfClassMacro(ImpactMetric)
itkStaticConstMacro(MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
unsigned long m_CurrentIteration
itkStaticConstMacro(FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
itk::SingleValuedCostFunction ITKBaseType
typename ElastixType::MovingImageType MovingImageType
typename ElastixType::RegistrationBaseType RegistrationType
A semantic similarity metric for multimodal image registration based on deep learning features.
std::vector< bool > GetBooleanVectorFromString(std::string valueStr, bool defaultValue)
Parse a string into typed values using a custom delimiter.
std::vector< std::vector< T > > GroupByDimensions(const std::vector< T > &values, const std::vector< unsigned int > &dimensions)
Groups a flat vector into sub-vectors based on given dimension sizes.


Generated on 1774142652 for elastix by doxygen 1.15.0 elastix logo