1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300
|
/*=========================================================================
*
* 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.
*
*=========================================================================*/
#ifndef itkTimeVaryingVelocityFieldIntegrationImageFilter_hxx
#define itkTimeVaryingVelocityFieldIntegrationImageFilter_hxx
#include "itkImageRegionIteratorWithIndex.h"
#include "itkVectorLinearInterpolateImageFunction.h"
namespace itk
{
/*
* TimeVaryingVelocityFieldIntegrationImageFilter class definitions
*/
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField,
TDisplacementField>::TimeVaryingVelocityFieldIntegrationImageFilter()
{
this->m_LowerTimeBound = 0.0, this->m_UpperTimeBound = 1.0, this->m_NumberOfIntegrationSteps = 100;
this->m_NumberOfTimePoints = 0;
this->SetNumberOfRequiredInputs(1);
if (InputImageDimension - 1 != OutputImageDimension)
{
itkExceptionMacro("The time-varying velocity field (input) should have "
<< "dimensionality of 1 greater than the deformation field (output). ");
}
using DefaultVelocityFieldInterpolatorType =
VectorLinearInterpolateImageFunction<TimeVaryingVelocityFieldType, ScalarType>;
typename DefaultVelocityFieldInterpolatorType::Pointer velocityFieldInterpolator =
DefaultVelocityFieldInterpolatorType::New();
this->m_VelocityFieldInterpolator = velocityFieldInterpolator;
using DefaultDisplacementFieldInterpolatorType =
VectorLinearInterpolateImageFunction<DisplacementFieldType, ScalarType>;
typename DefaultDisplacementFieldInterpolatorType::Pointer deformationFieldInterpolator =
DefaultDisplacementFieldInterpolatorType::New();
this->m_DisplacementFieldInterpolator = deformationFieldInterpolator;
this->DynamicMultiThreadingOn();
}
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
void
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField,
TDisplacementField>::GenerateOutputInformation()
{
const TimeVaryingVelocityFieldType * input = this->GetInput();
DisplacementFieldType * output = this->GetOutput();
this->m_NumberOfTimePoints = input->GetLargestPossibleRegion().GetSize()[OutputImageDimension];
if (!input || !output)
{
return;
}
//
// The ImageBase::CopyInformation() method ca not be used here
// because these two images have different dimensions. Therefore
// the individual elements must be copied for the common dimensions.
//
using SizeType = typename DisplacementFieldType::SizeType;
using SpacingType = typename DisplacementFieldType::SpacingType;
using OriginType = typename DisplacementFieldType::PointType;
using DirectionType = typename DisplacementFieldType::DirectionType;
SizeType size;
SpacingType spacing;
OriginType origin;
DirectionType direction;
using InputSizeType = typename TimeVaryingVelocityFieldType::SizeType;
using InputSpacingType = typename TimeVaryingVelocityFieldType::SpacingType;
using InputOriginType = typename TimeVaryingVelocityFieldType::PointType;
using InputDirectionType = typename TimeVaryingVelocityFieldType::DirectionType;
using InputRegionType = typename TimeVaryingVelocityFieldType::RegionType;
const InputSpacingType & inputSpacing = input->GetSpacing();
const InputOriginType & inputOrigin = input->GetOrigin();
const InputDirectionType & inputDirection = input->GetDirection();
const InputRegionType requestedRegion = input->GetRequestedRegion();
const InputSizeType requestedSize = requestedRegion.GetSize();
for (unsigned int i = 0; i < OutputImageDimension; ++i)
{
size[i] = requestedSize[i];
spacing[i] = inputSpacing[i];
origin[i] = inputOrigin[i];
for (unsigned int j = 0; j < OutputImageDimension; ++j)
{
direction[i][j] = inputDirection[i][j];
}
}
output->SetOrigin(origin);
output->SetSpacing(spacing);
output->SetDirection(direction);
output->SetRegions(size);
}
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
void
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField,
TDisplacementField>::BeforeThreadedGenerateData()
{
this->m_VelocityFieldInterpolator->SetInputImage(this->GetInput());
this->m_NumberOfTimePoints = this->GetInput()->GetLargestPossibleRegion().GetSize()[InputImageDimension - 1];
if (!this->m_InitialDiffeomorphism.IsNull())
{
this->m_DisplacementFieldInterpolator->SetInputImage(this->m_InitialDiffeomorphism);
}
}
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
void
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField, TDisplacementField>::
DynamicThreadedGenerateData(const OutputRegionType & region)
{
if (Math::ExactlyEquals(this->m_LowerTimeBound, this->m_UpperTimeBound) || this->m_NumberOfIntegrationSteps == 0)
{
this->GetOutput()->FillBuffer(typename DisplacementFieldType::PixelType{});
return;
}
const TimeVaryingVelocityFieldType * inputField = this->GetInput();
typename DisplacementFieldType::Pointer outputField = this->GetOutput();
for (ImageRegionIteratorWithIndex<DisplacementFieldType> It(outputField, region); !It.IsAtEnd(); ++It)
{
PointType point;
outputField->TransformIndexToPhysicalPoint(It.GetIndex(), point);
VectorType displacement = this->IntegrateVelocityAtPoint(point, inputField);
It.Set(displacement);
}
}
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
auto
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField, TDisplacementField>::IntegrateVelocityAtPoint(
const PointType & initialSpatialPoint,
const TimeVaryingVelocityFieldType * inputField) -> VectorType
{
// Solve the initial value problem using fourth-order Runge-Kutta
// y' = f(t, y), y(t_0) = y_0
VectorType zeroVector;
zeroVector.Fill(0.0);
// Initial conditions
VectorType displacement = zeroVector;
if (!this->m_InitialDiffeomorphism.IsNull())
{
if (this->m_DisplacementFieldInterpolator->IsInsideBuffer(initialSpatialPoint))
{
displacement = this->m_DisplacementFieldInterpolator->Evaluate(initialSpatialPoint);
}
}
// Perform the integration.
// With TimeBoundsAsRates On, we need to map the time dimension of the input image to the
// normalized domain of [0,1].
RealType timeOrigin = 0.;
RealType timeScale = 1.;
if (m_TimeBoundsAsRates)
{
typename TimeVaryingVelocityFieldType::PointType spaceTimeOrigin = inputField->GetOrigin();
using RegionType = typename TimeVaryingVelocityFieldType::RegionType;
RegionType region = inputField->GetLargestPossibleRegion();
typename RegionType::IndexType lastIndex = region.GetIndex();
typename RegionType::SizeType size = region.GetSize();
for (unsigned int d = 0; d < InputImageDimension; ++d)
{
lastIndex[d] += (size[d] - 1);
}
typename TimeVaryingVelocityFieldType::PointType spaceTimeEnd;
inputField->TransformIndexToPhysicalPoint(lastIndex, spaceTimeEnd);
timeOrigin = spaceTimeOrigin[InputImageDimension - 1];
const RealType timeEnd = spaceTimeEnd[InputImageDimension - 1];
timeScale = timeEnd - timeOrigin;
}
RealType timePointInImage = timeOrigin + this->m_LowerTimeBound * timeScale;
// Calculate the delta time used for integration
const RealType deltaTime =
(this->m_UpperTimeBound - this->m_LowerTimeBound) / static_cast<RealType>(this->m_NumberOfIntegrationSteps);
const RealType deltaTimeInImage = timeScale * deltaTime;
for (unsigned int n = 0; n < this->m_NumberOfIntegrationSteps; ++n)
{
typename TimeVaryingVelocityFieldType::PointType x1;
typename TimeVaryingVelocityFieldType::PointType x2;
typename TimeVaryingVelocityFieldType::PointType x3;
typename TimeVaryingVelocityFieldType::PointType x4;
for (unsigned int d = 0; d < OutputImageDimension; ++d)
{
x1[d] = initialSpatialPoint[d] + displacement[d];
x2[d] = x1[d];
x3[d] = x1[d];
x4[d] = x1[d];
}
x1[OutputImageDimension] = timePointInImage;
x2[OutputImageDimension] = timePointInImage + 0.5 * deltaTimeInImage;
x3[OutputImageDimension] = timePointInImage + 0.5 * deltaTimeInImage;
x4[OutputImageDimension] = timePointInImage + deltaTimeInImage;
VectorType f1 = zeroVector;
if (this->m_VelocityFieldInterpolator->IsInsideBuffer(x1))
{
f1 = this->m_VelocityFieldInterpolator->Evaluate(x1);
for (unsigned int jj = 0; jj < OutputImageDimension; ++jj)
{
x2[jj] += f1[jj] * deltaTime * 0.5;
}
}
VectorType f2 = zeroVector;
if (this->m_VelocityFieldInterpolator->IsInsideBuffer(x2))
{
f2 = this->m_VelocityFieldInterpolator->Evaluate(x2);
for (unsigned int jj = 0; jj < OutputImageDimension; ++jj)
{
x3[jj] += f2[jj] * deltaTime * 0.5;
}
}
VectorType f3 = zeroVector;
if (this->m_VelocityFieldInterpolator->IsInsideBuffer(x3))
{
f3 = this->m_VelocityFieldInterpolator->Evaluate(x3);
for (unsigned int jj = 0; jj < OutputImageDimension; ++jj)
{
x4[jj] += f3[jj] * deltaTime;
}
}
VectorType f4 = zeroVector;
if (this->m_VelocityFieldInterpolator->IsInsideBuffer(x4))
{
f4 = this->m_VelocityFieldInterpolator->Evaluate(x4);
}
for (unsigned int jj = 0; jj < OutputImageDimension; ++jj)
{
x1[jj] += deltaTime / 6.0 * (f1[jj] + 2.0 * f2[jj] + 2.0 * f3[jj] + f4[jj]);
displacement[jj] = x1[jj] - initialSpatialPoint[jj];
}
timePointInImage += deltaTimeInImage;
}
return displacement;
}
template <typename TTimeVaryingVelocityField, typename TDisplacementField>
void
TimeVaryingVelocityFieldIntegrationImageFilter<TTimeVaryingVelocityField, TDisplacementField>::PrintSelf(
std::ostream & os,
Indent indent) const
{
Superclass::PrintSelf(os, indent);
os << indent << "VelocityFieldInterpolator: " << this->m_VelocityFieldInterpolator << std::endl;
os << indent << "LowerTimeBound: " << this->m_LowerTimeBound << std::endl;
os << indent << "UpperTimeBound: " << this->m_UpperTimeBound << std::endl;
os << indent << "NumberOfIntegrationSteps: " << this->m_NumberOfIntegrationSteps << std::endl;
itkPrintSelfObjectMacro(InitialDiffeomorphism);
itkPrintSelfObjectMacro(DisplacementFieldInterpolator);
}
} // end namespace itk
#endif
|