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
|
/*=========================================================================
Program: VMTK
Module: $RCSfile: vtkvmtkVesselEnhancingDiffusionImageFilter.h,v $
Language: C++
Date: $Date: 2006/04/06 16:48:25 $
Version: $Revision: 1.2 $
Copyright (c) Luca Antiga, David Steinman. All rights reserved.
See LICENCE file for details.
Portions of this code are covered under the VTK copyright.
See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm
for details.
Portions of this code are covered under the ITK copyright.
See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm
for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
// .NAME vtkvmtkVesselEnhancingDiffusionImageFilter - Wrapper class around itk::VesselEnhancingDiffusionImageFilter
// .SECTION Description
// vtkvmtkVesselEnhancingDiffusionImageFilter
#ifndef __vtkvmtkVesselEnhancingDiffusionImageFilter_h
#define __vtkvmtkVesselEnhancingDiffusionImageFilter_h
#include "vtkvmtkITKImageToImageFilterFF.h"
#include "itkAnisotropicDiffusionVesselEnhancementImageFilter.h"
#include "itkMultiScaleHessianBasedMeasureImageFilter.h"
#include "itkHessianSmoothed3DToVesselnessMeasureImageFilter.h"
#include "vtkvmtkWin32Header.h"
class VTK_VMTK_SEGMENTATION_EXPORT vtkvmtkVesselEnhancingDiffusionImageFilter : public vtkvmtkITKImageToImageFilterFF
{
public:
static vtkvmtkVesselEnhancingDiffusionImageFilter *New();
vtkTypeRevisionMacro(vtkvmtkVesselEnhancingDiffusionImageFilter, vtkvmtkITKImageToImageFilterFF);
void SetSigmaMin(double value)
{
DelegateITKInputMacro(GetMultiScaleVesselnessFilter()->SetSigmaMinimum,value);
}
double GetSigmaMin()
{
DelegateITKOutputMacro(GetMultiScaleVesselnessFilter()->GetSigmaMinimum);
}
void SetSigmaMax(double value)
{
DelegateITKInputMacro(GetMultiScaleVesselnessFilter()->SetSigmaMaximum,value);
}
double GetSigmaMax()
{
DelegateITKOutputMacro(GetMultiScaleVesselnessFilter()->GetSigmaMaximum);
}
void SetNumberOfSigmaSteps(int value)
{
DelegateITKInputMacro(GetMultiScaleVesselnessFilter()->SetNumberOfSigmaSteps,value);
}
int GetNumberOfSigmaSteps()
{
DelegateITKOutputMacro(GetMultiScaleVesselnessFilter()->GetNumberOfSigmaSteps);
}
void SetSigmaStepMethodToEquispaced()
{
this->GetImageFilterPointer()->GetMultiScaleVesselnessFilter()->SetSigmaStepMethodToEquispaced();
this->Modified();
}
void SetSigmaStepMethodToLogarithmic()
{
this->GetImageFilterPointer()->GetMultiScaleVesselnessFilter()->SetSigmaStepMethodToLogarithmic();
this->Modified();
}
void SetTimeStep(double value)
{
DelegateITKInputMacro(SetTimeStep,value);
}
double GetTimeStep()
{
DelegateITKOutputMacro(GetTimeStep);
}
void SetEpsilon(double value)
{
DelegateITKInputMacro(SetEpsilon,value);
}
double GetEpsilon()
{
DelegateITKOutputMacro(GetEpsilon);
}
void SetWStrength(double value)
{
DelegateITKInputMacro(SetWStrength,value);
}
double GetWStrength()
{
DelegateITKOutputMacro(GetWStrength);
}
void SetSensitivity(double value)
{
DelegateITKInputMacro(SetSensitivity,value);
}
double GetSensitivity()
{
DelegateITKOutputMacro(GetSensitivity);
}
void SetNumberOfIterations(int value)
{
DelegateITKInputMacro(SetNumberOfIterations,value);
}
int GetNumberOfIterations()
{
DelegateITKOutputMacro(GetNumberOfIterations);
}
void SetAlpha(double value)
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->SetAlpha(value);
this->Modified();
}
}
double GetAlpha()
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
return dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->GetAlpha();
}
else
{
vtkErrorMacro ( << this->GetClassName() << " Error getting method. Dynamic cast returned 0" );
return 0.0;
}
}
void SetBeta(double value)
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->SetBeta(value);
this->Modified();
}
}
double GetBeta()
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
return dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->GetBeta();
}
else
{
vtkErrorMacro ( << this->GetClassName() << " Error getting method. Dynamic cast returned 0" );
return 0.0;
}
}
void SetGamma(double value)
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->SetGamma(value);
this->Modified();
}
}
double GetGamma()
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
return dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->GetGamma();
}
else
{
vtkErrorMacro ( << this->GetClassName() << " Error getting method. Dynamic cast returned 0" );
return 0.0;
}
}
void SetC(double value)
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->SetC(value);
this->Modified();
}
}
double GetC()
{
ImageFilterType* tempFilter = dynamic_cast<ImageFilterType*>(this->m_Filter.GetPointer());
if (tempFilter)
{
return dynamic_cast<VesselnessFilterType*>(tempFilter->GetMultiScaleVesselnessFilter()->GetHessianToMeasureFilter())->GetC();
}
else
{
vtkErrorMacro ( << this->GetClassName() << " Error getting method. Dynamic cast returned 0" );
return 0.0;
}
}
void SetNumberOfDiffusionSubIterations(int value)
{
DelegateITKInputMacro(SetNumberOfDiffusionSubIterations,value);
}
int GetNumberOfDiffusionSubIterations()
{
DelegateITKOutputMacro(GetNumberOfDiffusionSubIterations);
}
protected:
//BTX
typedef itk::AnisotropicDiffusionVesselEnhancementImageFilter<Superclass::InputImageType, Superclass::OutputImageType> ImageFilterType;
typedef ImageFilterType::VesselnessFilterType VesselnessFilterType;
vtkvmtkVesselEnhancingDiffusionImageFilter() : Superclass(ImageFilterType::New())
{
ImageFilterType* imageFilter = this->GetImageFilterPointer();
imageFilter->GetMultiScaleVesselnessFilter()->SetSigmaStepMethodToEquispaced();
}
~vtkvmtkVesselEnhancingDiffusionImageFilter() {};
ImageFilterType* GetImageFilterPointer() { return dynamic_cast<ImageFilterType*>(m_Filter.GetPointer()); }
//ETX
private:
vtkvmtkVesselEnhancingDiffusionImageFilter(const vtkvmtkVesselEnhancingDiffusionImageFilter&); // Not implemented.
void operator=(const vtkvmtkVesselEnhancingDiffusionImageFilter&); // Not implemented.
};
#endif
|