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
|
/*
//
// Copyright 2016 Google, Inc.
//
// Copyright 1997-2009 Torsten Rohlfing
//
// Copyright 2004-2011 SRI International
//
// This file is part of the Computational Morphometry Toolkit.
//
// http://www.nitrc.org/projects/cmtk/
//
// The Computational Morphometry Toolkit is free software: you can
// redistribute it and/or modify it under the terms of the GNU General Public
// License as published by the Free Software Foundation, either version 3 of
// the License, or (at your option) any later version.
//
// The Computational Morphometry Toolkit is distributed in the hope that it
// will be useful, but WITHOUT ANY WARRANTY; without even the implied
// warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU General Public License along
// with the Computational Morphometry Toolkit. If not, see
// <http://www.gnu.org/licenses/>.
//
// $Revision: 5436 $
//
// $LastChangedDate: 2018-12-10 19:01:20 -0800 (Mon, 10 Dec 2018) $
//
// $LastChangedBy: torstenrohlfing $
//
*/
#ifndef __cmtkFilterVolume_h_included_
#define __cmtkFilterVolume_h_included_
#include <list>
#include <cmtkconfig.h>
#include <Base/cmtkUniformVolume.h>
#include <Base/cmtkTypedArray.h>
#include <Base/cmtkUnits.h>
namespace
cmtk
{
/** \addtogroup Base */
//@{
/// Class for filtering volume images.
class FilterVolume
{
public:
/** Apply Gaussian filter.
*\param volume Input 3D image.
*\param width Width (standard deviation) of the Gaussian kernel.
*\param radius Filter radius in multiples of the filter width. Outside the
* radius (Euclidean distance) the filter is truncated.
*\param maskData Optional binary mask data array.
*\return A newly allocated TypedArray object that can be used, for
* example, to replace the one held by the input image. The data type of the
* array is identical to the input array.
*/
static TypedArray::SmartPtr GaussianFilter( const UniformVolume* volume, const Units::GaussianSigma& width, const Types::Coordinate radius = 1.0, const TypedArray* maskData = NULL );
/** Apply Coupe denoising filter.
*\param volume Input 3D image.
*\param beta Smoothing adjustment parameter
*\param windowRadius Distance from center voxel to outer edge of window
*\return A newly allocated TypedArray object that can be used, for
* example, to replace the one held by the input image. The data type of the
* array is identical to the input array.
*/
static TypedArray::SmartPtr CoupeFilter
( const UniformVolume* volume,
const int windowRadius,
const float beta = 0.5 );
/** Apply Torsten Rohlfing's single-image intensity-consistent Gaussian filter.
*\param volume Input 3D image.
*\param subjectData Pixel array of the individual grey image from this
* subject.
*\param maskData Optional binary mask data array.
*\param iFilterSigma Width (standard deviation of the Gaussian kernel.
*\param filterWidth Width (standard deviation of the Gaussian kernel.
*\param filterRadius Filter radius in multiples of the filter width.
* Outside the radius (Euclidean distance) the filter is truncated.
*\return A newly allocated TypedArray object that can be used, for
* example, to replace the one held by the input image. The data type of the
* array is identical to the input array.
*/
static TypedArray::SmartPtr RohlfingFilter
( const UniformVolume* volume, const TypedArray* subjectData,
const TypedArray* maskData, const Units::GaussianSigma& iFilterSigma,
const Units::GaussianSigma& filterWidth, const Types::Coordinate filterRadius );
/** Apply Colin Studholme's Gaussian filter with registration-based weights.
*\param volume Input 3D image.
*\param subjectData Pixel array of the individual grey image from this
* subject.
*\param averageData Pixel array of the population average grey image.
*\param maskData Optional binary mask data array.
*\param imgList List of pixel arrays from matched 3D images. The consistency
* between these pixel arrays determines the relative weights of the
* otherwise Gaussian kernel. Consult Colin's NeuroImage (2003) paper for
* details.
*\param binWidth Bin width of the intensity histogram used to quantify the
* matching between all individuals.
*\param filterWidth Width (standard deviation of the Gaussian kernel.
*\param filterRadius Filter radius in multiples of the filter width.
* Outside the radius (Euclidean distance) the filter is truncated.
*\return A newly allocated TypedArray object that can be used, for
* example, to replace the one held by the input image. The data type of the
* array is identical to the input array.
*/
static TypedArray::SmartPtr StudholmeFilter
( const UniformVolume* volume, const TypedArray* subjectData,
const TypedArray* averageData, const TypedArray* maskData,
std::list<TypedArray::SmartPtr> imgList, const Types::DataItem binWidth,
const Units::GaussianSigma& filterWidth, const Types::Coordinate filterRadius );
/** Apply Colin Studholme's Gaussian filter using multiple time points.
*\param volume Input 3D image.
*\param subjectData List of pixel arrays of the individual grey images from
* this subject. Each list item corresponds to one time point. Ideally, there
* are two items in the list, which correspond to the time points between
* which the Jacobian map in "volume" was computed.
*\param averageData Pixel array of the population average grey image.
*\param maskData Optional binary mask data array.
*\param imgList List of pixel arrays from matched 3D images. The consistency
* between these pixel arrays determines the relative weights of the
* otherwise Gaussian kernel. Consult Colin's NeuroImage (2003) paper for
* details.
*\param binWidth Bin width of the intensity histogram used to quantify the
* matching between all individuals.
*\param filterWidth Width (standard deviation of the Gaussian kernel.
*\param filterRadius Filter radius in multiples of the filter width.
* Outside the radius (Euclidean distance) the filter is truncated.
*\return A newly allocated TypedArray object that can be used, for
* example, to replace the one held by the input image. The data type of the
* array is identical to the input array.
*/
static TypedArray::SmartPtr StudholmeFilter
( const UniformVolume* volume,
std::list<TypedArray::SmartPtr> subjectData,
const TypedArray* averageData, const TypedArray* maskData,
std::list<TypedArray::SmartPtr> imgList, const Types::DataItem binWidth,
const Units::GaussianSigma& filterWidth, const Types::Coordinate filterRadius );
};
//@}
} // namespace cmtk
#endif // #ifndef __cmtkFilterVolume_h_included_
|