Skip to content

Commit 6a5060b

Browse files
blowekamphjmjohnson
authored andcommitted
ENH: Restore integer support for laplacian sharpening
Revert to using the a laplacian neighborhood operator over the recursive laplacian filter which requires real pixel types for input and output.
1 parent 14816ad commit 6a5060b

File tree

3 files changed

+37
-10
lines changed

3 files changed

+37
-10
lines changed

Examples/Filtering/LaplacianSharpeningImageFilter.cxx

Lines changed: 3 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -34,21 +34,19 @@ main(int argc, char * argv[])
3434
const char * inputFilename = argv[1];
3535
const char * outputFilename = argv[2];
3636

37-
using PixelType = float;
3837
using CharPixelType = unsigned char;
3938
constexpr unsigned int Dimension = 2;
4039

41-
using ImageType = itk::Image<PixelType, Dimension>;
4240
using CharImageType = itk::Image<CharPixelType, Dimension>;
4341

44-
using ReaderType = itk::ImageFileReader<ImageType>;
42+
using ReaderType = itk::ImageFileReader<CharImageType>;
4543
using WriterType = itk::ImageFileWriter<CharImageType>;
4644

4745
using RescaleFilter =
48-
itk::RescaleIntensityImageFilter<ImageType, CharImageType>;
46+
itk::RescaleIntensityImageFilter<CharImageType, CharImageType>;
4947

5048
using LaplacianSharpeningFilter =
51-
itk::LaplacianSharpeningImageFilter<ImageType, ImageType>;
49+
itk::LaplacianSharpeningImageFilter<CharImageType, CharImageType>;
5250

5351

5452
// Setting the IO

Modules/Filtering/ImageFeature/include/itkLaplacianSharpeningImageFilter.h

Lines changed: 0 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -93,8 +93,6 @@ class ITK_TEMPLATE_EXPORT LaplacianSharpeningImageFilter : public ImageToImageFi
9393
#ifdef ITK_USE_CONCEPT_CHECKING
9494
// Begin concept checking
9595
itkConceptMacro(SameDimensionCheck, (Concept::SameDimension<InputImageDimension, ImageDimension>));
96-
itkConceptMacro(InputPixelTypeIsFloatingPointCheck, (Concept::IsFloatingPoint<InputPixelType>));
97-
itkConceptMacro(OutputPixelTypeIsFloatingPointCheck, (Concept::IsFloatingPoint<OutputPixelType>));
9896
// End concept checking
9997
#endif
10098

Modules/Filtering/ImageFeature/include/itkLaplacianSharpeningImageFilter.hxx

Lines changed: 34 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -69,10 +69,41 @@ LaplacianSharpeningImageFilter<TInputImage, TOutputImage>::GenerateData()
6969

7070
using RealImageType = Image<RealType, ImageDimension>;
7171

72-
using LaplacianImageFilter = LaplacianImageFilter<InputImageType, RealImageType>;
73-
typename LaplacianImageFilter::Pointer laplacianFilter = LaplacianImageFilter::New();
72+
// Create the Laplacian operator
73+
LaplacianOperator<RealType, ImageDimension> oper;
74+
double s[ImageDimension];
75+
for (unsigned int i = 0; i < ImageDimension; ++i)
76+
{
77+
if (localInput->GetSpacing()[i] == 0.0)
78+
{
79+
itkExceptionMacro("Image spacing cannot be zero");
80+
}
81+
else if (this->m_UseImageSpacing)
82+
{
83+
s[i] = 1.0 / localInput->GetSpacing()[i];
84+
}
85+
else
86+
{
87+
s[i] = 1.0;
88+
}
89+
}
90+
oper.SetDerivativeScalings(s);
91+
oper.CreateOperator();
92+
// Calculate the Laplacian filtered image
93+
94+
// do calculations in floating point
95+
using RealImageType = Image<RealType, ImageDimension>;
96+
using NOIF = NeighborhoodOperatorImageFilter<InputImageType, RealImageType>;
97+
ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
98+
99+
auto laplacianFilter = NOIF::New();
100+
laplacianFilter->OverrideBoundaryCondition(static_cast<typename NOIF::ImageBoundaryConditionPointerType>(&nbc));
101+
102+
//
103+
// set up the mini-pipeline
104+
//
105+
laplacianFilter->SetOperator(oper);
74106
laplacianFilter->SetInput(localInput);
75-
laplacianFilter->SetUseImageSpacing(m_UseImageSpacing);
76107
laplacianFilter->Update();
77108

78109
auto filteredMinMaxFilter = MinimumMaximumImageFilter<RealImageType>::New();

0 commit comments

Comments
 (0)