Skip to content

Commit

Permalink
ENH: Import itkBSplineGradientImageFilter
Browse files Browse the repository at this point in the history
  • Loading branch information
thewtex committed Jul 30, 2018
1 parent 31dbeee commit 7fe2172
Show file tree
Hide file tree
Showing 5 changed files with 315 additions and 2 deletions.
4 changes: 2 additions & 2 deletions include/itkBSplineApproximationGradientImageFilter.h
Original file line number Diff line number Diff line change
Expand Up @@ -43,8 +43,8 @@ namespace itk
*
* \ingroup BSplineGradient
*/
template< class TInputImage, class TOutputValueType >
class BSplineApproximationGradientImageFilter:
template< typename TInputImage, typename TOutputValueType >
class ITK_TEMPLATE_EXPORT BSplineApproximationGradientImageFilter:
public ImageToImageFilter< TInputImage,
Image< CovariantVector< TOutputValueType,
TInputImage::ImageDimension >,
Expand Down
122 changes: 122 additions & 0 deletions include/itkBSplineGradientImageFilter.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,122 @@
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* 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
*
* http://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 itkBSplineGradientImageFilter_h
#define itkBSplineGradientImageFilter_h

#include "itkBSplineInterpolateImageFunction.h"
#include "itkImageToImageFilter.h"
#include "itkImage.h"
#include "itkCovariantVector.h"

namespace itk
{
/** \class BSplineGradientImageFilter
*
* \brief Calculate image gradient from an interpolating B-spline fit.
*
* \ingroup GradientFilters
* \ingroup BSplineGradient
*
* \sa ApproximationBSplineGradientImageFilter
*/
template< typename TInputImage,
typename TOutputValueType = float,
typename TCoordRep = double,
typename TCoefficientType = double >
class ITK_TEMPLATE_EXPORT BSplineGradientImageFilter :
public ImageToImageFilter< TInputImage,
Image< CovariantVector< TOutputValueType,
TInputImage::ImageDimension >,
TInputImage::ImageDimension > >
{
public:
ITK_DISALLOW_COPY_AND_ASSIGN(BSplineGradientImageFilter);

/** Extract dimension from input image. */
itkStaticConstMacro(ImageDimension, unsigned int, TInputImage::ImageDimension);

/** Standard class typedefs. */
typedef BSplineGradientImageFilter Self;

/** Convenient typedefs for simplifying declarations. */
typedef TInputImage InputImageType;
typedef typename InputImageType::Pointer InputImagePointer;
typedef Image< CovariantVector<
TOutputValueType, ImageDimension >,
ImageDimension >
OutputImageType;

/** Standard class typedefs. */
typedef ImageToImageFilter< InputImageType, OutputImageType > Superclass;
typedef SmartPointer< Self > Pointer;
typedef SmartPointer< const Self > ConstPointer;

/** Method for creation through the object factory. */
itkNewMacro(Self);

/** Run-time type information (and related methods). */
itkTypeMacro(BSplineGradientImageFilter, ImageToImageFilter);

/** Image typedef support. */
typedef typename InputImageType::PixelType InputPixelType;
typedef TOutputValueType OutputValueType;
typedef CovariantVector<
OutputValueType, itkGetStaticConstMacro(OutputImageDimension) >
OutputPixelType;
typedef typename OutputImageType::RegionType OutputImageRegionType;

#ifdef ITK_USE_CONCEPT_CHECKING
/** Begin concept checking */
itkConceptMacro( InputConvertibleToOutputCheck,
( Concept::Convertible< InputPixelType, OutputValueType > ) );
itkConceptMacro( OutputHasNumericTraitsCheck,
( Concept::HasNumericTraits< OutputValueType > ) );
/** End concept checking */
#endif

/** Typedefs for the interpolator. */
typedef TCoordRep CoordRepType;
typedef TCoefficientType CoefficientType;
typedef itk::BSplineInterpolateImageFunction< InputImageType, CoordRepType, CoefficientType >
InterpolatorType;
typedef typename InterpolatorType::Pointer InterpolatorPointerType;

protected:
BSplineGradientImageFilter();
virtual ~BSplineGradientImageFilter() {}

/** This filter requires all of the input image for now. Future
* implementation may be capable of streaming. */
void GenerateInputRequestedRegion() override;

/** This filter must produce all of its output at once. */
void EnlargeOutputRequestedRegion(DataObject *output) override;

void DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;

private:

};

} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#include "itkBSplineGradientImageFilter.hxx"
#endif

#endif
114 changes: 114 additions & 0 deletions include/itkBSplineGradientImageFilter.hxx
Original file line number Diff line number Diff line change
@@ -0,0 +1,114 @@
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* 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
*
* http://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 itkBSplineGradientImageFilter_hxx
#define itkBSplineGradientImageFilter_hxx

#include "itkBSplineGradientImageFilter.h"

#include "itkImageRegionIteratorWithIndex.h"

namespace itk
{

template< typename TInputImage, typename TOutputValueType,
typename TCoordRep, typename TCoefficientType >
BSplineGradientImageFilter< TInputImage, TOutputValueType,
TCoordRep, TCoefficientType >
::BSplineGradientImageFilter()
{
}


template< typename TInputImage, typename TOutputValueType,
typename TCoordRep, typename TCoefficientType >
void
BSplineGradientImageFilter< TInputImage, TOutputValueType,
TCoordRep, TCoefficientType >
::GenerateInputRequestedRegion()
{
// this filter requires the all of the input image to be in
// the buffer
InputImagePointer inputPtr = const_cast< InputImageType * >( this->GetInput() );

if ( inputPtr )
{
inputPtr->SetRequestedRegionToLargestPossibleRegion();
}
}


template< typename TInputImage, typename TOutputValueType,
typename TCoordRep, typename TCoefficientType >
void
BSplineGradientImageFilter< TInputImage, TOutputValueType,
TCoordRep, TCoefficientType >
::EnlargeOutputRequestedRegion( DataObject *output)
{
// this filter requires the all of the output image to be in
// the buffer
OutputImageType *imgData;

imgData = dynamic_cast< OutputImageType * >( output );
if ( imgData )
{
imgData->SetRequestedRegionToLargestPossibleRegion();
}
}


template< typename TInputImage, typename TOutputValueType,
typename TCoordRep, typename TCoefficientType >
void
BSplineGradientImageFilter< TInputImage, TOutputValueType,
TCoordRep, TCoefficientType >
::DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread)
{
typename OutputImageType::Pointer outputPtr = this->GetOutput();

if( !outputPtr )
{
return;
}

// this filter requires the all of the input image to be in
// the buffer
InputImagePointer inputPtr = const_cast< TInputImage * >( this->GetInput() );
InterpolatorPointerType interpolator = InterpolatorType::New();
// This will calculate the coefficients.
interpolator->SetInputImage( inputPtr );

typedef typename itk::ImageRegionIteratorWithIndex< OutputImageType > IteratorType;
typename OutputImageType::IndexType index;
typename InterpolatorType::ContinuousIndexType contIndex;
unsigned int i;
IteratorType it( outputPtr, outputRegionForThread );

for( it.GoToBegin(); !it.IsAtEnd(); ++it )
{
index = it.GetIndex();
for( i = 0; i < ImageDimension; ++i )
{
contIndex[i] = static_cast< CoordRepType >( index[i] );
}
it.Set( interpolator->EvaluateDerivativeAtContinuousIndex( contIndex ) );
}
}

} // end namespace itk

#endif
8 changes: 8 additions & 0 deletions test/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@ set(BSplineGradientTests
itkImageToImageOfVectorsFilterTest.cxx
itkBSplineScatteredDataPointSetToGradientImageFilterTest.cxx
itkBSplineApproximationGradientImageFilterTest.cxx
itkBSplineGradientImageFilterTest.cxx
)

CreateTestDriver(BSplineGradient "${BSplineGradient-Test_LIBRARIES}" "${BSplineGradientTests}")
Expand Down Expand Up @@ -36,3 +37,10 @@ itk_add_test(NAME itkBSplineApproximationGradientImageFilterTest
DATA{Input/foot.mha}
${ITK_TEST_OUTPUT_DIR}/itkBSplineApproximationGradientImageFilterTest
)

itk_add_test(NAME itkBSplineGradientImageFilterTest
COMMAND BSplineGradientTestDriver
itkBSplineGradientImageFilterTest
DATA{Input/foot.mha}
${ITK_TEST_OUTPUT_DIR}/itkBSplineGradientImageFilterTestGradientMagnitude.mha
)
69 changes: 69 additions & 0 deletions test/itkBSplineGradientImageFilterTest.cxx
Original file line number Diff line number Diff line change
@@ -0,0 +1,69 @@
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* 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
*
* http://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.
*
*=========================================================================*/
#include "itkVectorMagnitudeImageFilter.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"

#include "itkBSplineGradientImageFilter.h"

int itkBSplineGradientImageFilterTest(int argc, char *argv[])
{
if ( argc < 3 )
{
std::cerr << "Usage: " << argv[0];
std::cerr << " inputImage outputImage ";
std::cerr << std::endl;
return EXIT_FAILURE;
}

const unsigned int Dimension = 2;
typedef float PixelType;
typedef itk::Image< PixelType, Dimension > InputImageType;

typedef itk::ImageFileReader< InputImageType > ReaderType;
ReaderType::Pointer reader = ReaderType::New();

typedef itk::BSplineGradientImageFilter< InputImageType, PixelType > FilterType;
FilterType::Pointer filter = FilterType::New();
filter->SetInput( reader->GetOutput() );
reader->SetFileName( argv[1] );

typedef FilterType::OutputImageType GradientImageType;
typedef itk::VectorMagnitudeImageFilter< GradientImageType, InputImageType >
GradientMagnitudeFilterType;
GradientMagnitudeFilterType::Pointer gradientMagnitude = GradientMagnitudeFilterType::New();
gradientMagnitude->SetInput( filter->GetOutput() );

typedef itk::ImageFileWriter< InputImageType > WriterType;
WriterType::Pointer writer = WriterType::New();
writer->SetInput( gradientMagnitude->GetOutput() );
writer->SetFileName( argv[2] );

try
{
writer->Update();
}
catch( itk::ExceptionObject& ex )
{
std::cerr << "Exception caught!" << std::endl;
std::cerr << ex << std::endl;
return EXIT_FAILURE;
}

return EXIT_SUCCESS;
}

0 comments on commit 7fe2172

Please sign in to comment.