/*========================================================================= * * Copyright NumFOCUS * * Licensed under the Apache License, Version 2.0 (the "License"); * you may not use this file except in compliance with the License. * You may obtain a copy of the License at * * https://www.apache.org/licenses/LICENSE-2.0.txt * * Unless required by applicable law or agreed to in writing, software * distributed under the License is distributed on an "AS IS" BASIS, * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. * See the License for the specific language governing permissions and * limitations under the License. * *=========================================================================*/ #define ITK_TEMPLATE_EXPLICIT_ObjectToObjectOptimizerBaseTemplate #include "itkObjectToObjectOptimizerBase.h" #include "itkMultiThreaderBase.h" namespace itk { ITK_GCC_PRAGMA_DIAG_PUSH() ITK_GCC_PRAGMA_DIAG(ignored "-Wattributes") template ObjectToObjectOptimizerBaseTemplate::ObjectToObjectOptimizerBaseTemplate() { this->m_Metric = nullptr; this->m_CurrentIteration = 0; this->m_NumberOfIterations = 100; this->m_CurrentMetricValue = 0; // Initialize, but w/out calling SetNumberOfWorkUnits, to avoid // valgrind warning. this->m_NumberOfWorkUnits = MultiThreaderBase::GetGlobalDefaultNumberOfThreads(); this->m_ScalesAreIdentity = false; this->m_WeightsAreIdentity = true; this->m_DoEstimateScales = true; } template ObjectToObjectOptimizerBaseTemplate::~ObjectToObjectOptimizerBaseTemplate() = default; template void ObjectToObjectOptimizerBaseTemplate::PrintSelf(std::ostream & os, Indent indent) const { Superclass::PrintSelf(os, indent); itkPrintSelfObjectMacro(Metric); os << indent << "NumberOfWorkUnits: " << static_cast::PrintType>(m_NumberOfWorkUnits) << std::endl; os << indent << "CurrentIteration: " << static_cast::PrintType>(m_CurrentIteration) << std::endl; os << indent << "NumberOfIterations: " << static_cast::PrintType>(m_NumberOfIterations) << std::endl; os << indent << "CurrentMetricValue: " << static_cast::PrintType>(m_CurrentMetricValue) << std::endl; os << indent << "Scales: " << static_cast::PrintType>(m_Scales) << std::endl; os << indent << "Weights: " << static_cast::PrintType>(m_Weights) << std::endl; os << indent << "ScalesAreIdentity: " << (m_ScalesAreIdentity ? "On" : "Off") << std::endl; itkPrintSelfObjectMacro(ScalesEstimator); os << indent << "WeightsAreIdentity: " << (m_WeightsAreIdentity ? "On" : "Off") << std::endl; os << indent << "DoEstimateScales: " << (m_DoEstimateScales ? "On" : "Off") << std::endl; } template void ObjectToObjectOptimizerBaseTemplate::SetNumberOfWorkUnits(ThreadIdType number) { if (number < 1) { itkExceptionMacro("Number of work units must be > 0"); } if (number != this->m_NumberOfWorkUnits) { this->m_NumberOfWorkUnits = number; this->Modified(); } } template void ObjectToObjectOptimizerBaseTemplate::StartOptimization( bool itkNotUsed(doOnlyInitialization)) { /* Validate some settings */ if (this->m_Metric.IsNull()) { itkExceptionMacro("m_Metric must be set."); } /* Estimate the parameter scales if requested. */ if (this->m_DoEstimateScales && this->m_ScalesEstimator.IsNotNull()) { ScalesType scales; this->m_ScalesEstimator->EstimateScales(scales); this->SetScales(scales); itkDebugMacro("Estimated scales = " << this->m_Scales); } /* Verify m_Scales. If m_Scales hasn't been set, initialize to all 1's. */ using SValueType = typename ScalesType::ValueType; if (this->GetScalesInitialized()) { if (this->m_Scales.Size() != this->m_Metric->GetNumberOfLocalParameters()) { itkExceptionMacro("Size of scales (" << this->m_Scales.Size() << ") must equal number of local parameters (" << this->m_Metric->GetNumberOfLocalParameters() << ")."); } /* Check that all values in m_Scales are > machine epsilon, to avoid * division by zero/epsilon. * Also check if scales are identity. */ using SizeType = typename ScalesType::size_type; this->m_ScalesAreIdentity = true; for (SizeType i = 0; i < this->m_Scales.Size(); ++i) { if (this->m_Scales[i] <= NumericTraits::epsilon()) { itkExceptionMacro("m_Scales values must be > epsilon." << this->m_Scales); } /* Check if the scales are identity. Consider to be identity if * within a tolerance, to allow for automatically estimated scales * that may not be exactly 1.0 when in priciniple they should be. */ SValueType difference = itk::Math::abs(NumericTraits::OneValue() - this->m_Scales[i]); auto tolerance = static_cast(0.01); if (difference > tolerance) { this->m_ScalesAreIdentity = false; break; } } } else { // Initialize scales to identity m_Scales.SetSize(this->m_Metric->GetNumberOfLocalParameters()); m_Scales.Fill(NumericTraits::OneValue()); this->m_ScalesAreIdentity = true; } /* Verify m_Weights. */ using SValueType = typename ScalesType::ValueType; if (this->m_Weights.Size() > 0) { if (this->m_Weights.Size() != this->m_Metric->GetNumberOfLocalParameters()) { itkExceptionMacro("Size of weights (" << this->m_Weights.Size() << ") must equal number of local parameters (" << this->m_Metric->GetNumberOfLocalParameters() << ")."); } /* Check if they are identity within tolerance. */ using SizeType = typename ScalesType::size_type; this->m_WeightsAreIdentity = true; for (SizeType i = 0; i < this->m_Weights.Size(); ++i) { SValueType difference = itk::Math::abs(NumericTraits::OneValue() - this->m_Weights[i]); auto tolerance = static_cast(1e-4); if (difference > tolerance) { this->m_WeightsAreIdentity = false; break; } } } else { // Set weights to identity. But leave the array empty. this->m_WeightsAreIdentity = true; } } template const typename ObjectToObjectOptimizerBaseTemplate::ParametersType & ObjectToObjectOptimizerBaseTemplate::GetCurrentPosition() const { if (this->m_Metric.IsNull()) { itkExceptionMacro("m_Metric has not been assigned. Cannot get parameters."); } return this->m_Metric->GetParameters(); } template const typename ObjectToObjectOptimizerBaseTemplate::MeasureType & ObjectToObjectOptimizerBaseTemplate::GetValue() const { return this->GetCurrentMetricValue(); } template bool ObjectToObjectOptimizerBaseTemplate::GetScalesInitialized() const { return m_Scales.Size() > 0; } template class ITKOptimizersv4_EXPORT ObjectToObjectOptimizerBaseTemplate; template class ITKOptimizersv4_EXPORT ObjectToObjectOptimizerBaseTemplate; ITK_GCC_PRAGMA_DIAG_POP() /** Print enum values */ std::ostream & operator<<(std::ostream & out, const ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer value) { return out << [value] { switch (value) { case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::MAXIMUM_NUMBER_OF_ITERATIONS: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::MAXIMUM_NUMBER_OF_" "ITERATIONS"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::COSTFUNCTION_ERROR: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::COSTFUNCTION_" "ERROR"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::UPDATE_PARAMETERS_ERROR: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::UPDATE_PARAMETERS_" "ERROR"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::STEP_TOO_SMALL: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::STEP_TOO_SMALL"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::CONVERGENCE_CHECKER_PASSED: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::CONVERGENCE_" "CHECKER_PASSED"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer:: GRADIENT_MAGNITUDE_TOLEARANCE: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::GRADIENT_" "MAGNITUDE_TOLEARANCE"; case ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::OTHER_ERROR: return "itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer::OTHER_ERROR"; default: return "INVALID VALUE FOR itk::ObjectToObjectOptimizerBaseTemplateEnums::StopConditionObjectToObjectOptimizer"; } }(); } } // end namespace itk