Medical Imaging Interaction Toolkit  2016.11.0
Medical Imaging Interaction Toolkit
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
itkMultiGaussianImageSource.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef __itkMultiGaussianImageSource_h
29 #define __itkMultiGaussianImageSource_h
30 
31 #include "itkImageSource.h"
32 #include "itkNumericTraits.h"
33 #include "itkImageRegionIteratorWithIndex.h"
34 #include "itkImageFileWriter.h"
35 #include <itkMapContainer.h>
36 
37 
38 namespace itk
39 {
193 template< typename TOutputImage >
194 class ITK_EXPORT MultiGaussianImageSource:public ImageSource< TOutputImage >
195 {
196 public:
199  typedef ImageSource< TOutputImage > Superclass;
202 
205 
207  typedef typename TOutputImage::RegionType OutputImageRegionType;
208 
210  itkNewMacro(Self);
211 
213  typedef typename TOutputImage::SizeType SizeType;
214  typedef typename TOutputImage::IndexType IndexType;
215  typedef typename TOutputImage::SpacingType SpacingType;
217  typedef typename SizeType::SizeValueType SizeValueType;
218  typedef SizeValueType SizeValueArrayType[TOutputImage::ImageDimension];
219  typedef typename TOutputImage::SpacingValueType SpacingValueType;
220  typedef SpacingValueType SpacingValueArrayType[TOutputImage::ImageDimension];
221  typedef typename TOutputImage::PointValueType PointValueType;
222  typedef PointValueType PointValueArrayType[TOutputImage::ImageDimension];
223  typedef typename itk::ImageRegion<3> ::SizeValueType SizeRegionType;
225  typedef double RadiusType;
227  typedef std::vector<double> VectorType;
229  typedef Vector<double, TOutputImage::ImageDimension> ItkVectorType;
231  typedef ImageRegionIteratorWithIndex<TOutputImage> IteratorType;
234 
235  typedef MapContainer<unsigned int, PointType> MapContainerPoints;
236  typedef MapContainer<unsigned int, double> MapContainerRadius;
237 
238 
240  itkSetMacro(Size, SizeType);
241  virtual void SetSize(SizeValueArrayType sizeArray);
242  virtual const SizeValueType * GetSize() const;
244  itkSetMacro(Spacing, SpacingType);
245  virtual void SetSpacing(SpacingValueArrayType spacingArray);
246  virtual const SpacingValueType * GetSpacing() const;
248  itkSetMacro(Origin, PointType);
249  virtual void SetOrigin(PointValueArrayType originArray);
250  virtual const PointValueType * GetOrigin() const;
252  virtual unsigned int GetNumberOfGaussians() const;
254  virtual void SetNumberOfGausssians( unsigned int );
256  virtual RadiusType GetRadius() const;
257  virtual void SetRadius( RadiusType radius );
259  virtual const OutputImagePixelType GetMaxMeanValue() const;
261  virtual const IndexType GetSphereMidpoint() const;
263  virtual double MultiGaussianFunctionValueAtPoint(double , double, double);
266  virtual void AddGaussian( VectorType centerX, VectorType centerY, VectorType centerZ, VectorType sigmaX, VectorType sigmaY, VectorType sigmaZ, VectorType altitude);
268  virtual void CalculateTheMidpointAndTheMeanValueWithOctree();
270  virtual void CalculateMaxAndMinInSphere();
272  virtual const IndexType GetMaxValueIndexInSphere() const;
274  virtual const OutputImagePixelType GetMaxValueInSphere() const;
276  virtual const IndexType GetMinValueIndexInSphere() const;
278  virtual const OutputImagePixelType GetMinValueInSphere() const;
280  virtual void SetRegionOfInterest(ItkVectorType, ItkVectorType);
282  virtual void WriteXMLToTestTheCuboidInsideTheSphere();
284  virtual void CalculateEdgesInSphere( PointType globalCoordinateMidpointCuboid, PointType globalCoordinateMidpointSphere, double cuboidRadius, int level);
286  virtual double MultiGaussianFunctionValueAtCuboid(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
288  virtual void InsertPoints( PointType globalCoordinateMidpointCuboid, double cuboidRadius);
290  virtual void GenerateCuboidSegmentationInSphere( PointType globalCoordinateMidpointSphere );
292  virtual double FunctionPhi(double value);
294  virtual unsigned int IntersectTheSphere( PointType globalCoordinateMidpointCuboid, PointType globalCoordinateMidpointSphere, double sideLength);
296  void SetNormalDistributionValues();
297 
298 
301  itkSetClampMacro( Min, OutputImagePixelType,
302  NumericTraits< OutputImagePixelType >::NonpositiveMin(),
305  bool IsInImage(IndexType index);
306 
308  itkGetConstMacro(Min, OutputImagePixelType);
309 
312  itkSetClampMacro( Max, OutputImagePixelType,
313  NumericTraits< OutputImagePixelType >::NonpositiveMin(),
315 
317  itkGetConstMacro(Max, OutputImagePixelType);
318 
319 protected:
322  void PrintSelf(std::ostream & os, Indent indent) const;
323 
324  virtual void GenerateData();
325  virtual void GenerateOutputInformation();
326 
327 private:
328  MultiGaussianImageSource(const MultiGaussianImageSource &); //purposely not implemented
329  void operator=(const MultiGaussianImageSource &); //purposely not implemented
330 
331  SizeType m_Size; //size of the output image
332  SpacingType m_Spacing; //spacing
333  PointType m_Origin; //origin
334  OutputImagePixelType m_MaxValueInSphere; //maximal value in the wanted sphere
335  IndexType m_MaxValueIndexInSphere; //index of the maximal value in the wanted sphere
336  OutputImagePixelType m_MinValueInSphere; //minimal value in the wanted sphere
337  IndexType m_MinValueIndexInSphere; //index of the minimal value in the wanted sphere
338  unsigned int m_NumberOfGaussians; //number of Gaussians
339  RadiusType m_Radius; //radius of the sphere
340  unsigned int m_RadiusStepNumber; //number of steps to traverse the sphere radius
341  OutputImagePixelType m_MeanValue; //mean value in the wanted sphere
342  OutputImagePixelType m_ValueAtMidpoint; //value at the midpoint of the wanted sphere
343  IndexType m_SphereMidpoint; //midpoint of the wanted sphere
344  VectorType m_SigmaX; //deviation in the x-axis
345  VectorType m_SigmaY; //deviation in the y-axis
346  VectorType m_SigmaZ; //deviation in the z-axis
347  VectorType m_CenterX; //x-coordinate of the mean value of Gaussians
348  VectorType m_CenterY; //y-coordinate of the mean value of Gaussians
349  VectorType m_CenterZ; //z-coordinate of the mean value of Gaussians
350  VectorType m_Altitude; //amplitude
351  ItkVectorType m_RegionOfInterestMax; //maximal values for the coordinates in the region of interest
352  ItkVectorType m_RegionOfInterestMin; //minimal values for the coordinates in the region of interest
353  typename TOutputImage::PixelType m_Min; //minimum possible value
354  typename TOutputImage::PixelType m_Max; //maximum possible value
355  PointType m_GlobalCoordinate; //physical coordiante of the sphere midpoint
356  bool m_WriteMPS; //1 = write a MPS File to visualise the cuboid midpoints of one approximation of the sphere
357  MapContainerPoints m_Midpoints; //the midpoints of the cuboids
358  MapContainerRadius m_RadiusCuboid; //the radius ( = 0.5 * side length) of the cuboids (in the same order as the midpoints in m_Midpoints)
359  double m_Volume; //the volume of the body, that approximize the sphere
360  double m_NormalDistValues [410];//normal distribution values
361  double m_meanValueTemp; //= m_Volume * meanValue in each sphere
362  // The following variables are deprecated, and provided here just for
363  // backward compatibility. It use is discouraged.
364  mutable PointValueArrayType m_OriginArray;
365  mutable SpacingValueArrayType m_SpacingArray;
366 };
367 
368 } // end namespace itk
369 
370 #ifndef ITK_MANUAL_INSTANTIATION
371 #include "itkMultiGaussianImageSource.hxx"
372 #endif
373 
374 #endif
mitk::Point3D PointType
itk::SmartPointer< Self > Pointer
TOutputImage::PointValueType PointValueType
#define Max(x, y)
Definition: AnnotationP.h:41
TOutputImage::SpacingValueType SpacingValueType
Generate an 3-dimensional multigaussian image.
ImageRegionIteratorWithIndex< TOutputImage > IteratorType
TOutputImage::RegionType OutputImageRegionType
MapContainer< unsigned int, PointType > MapContainerPoints
MapContainer< unsigned int, double > MapContainerRadius
static T max(T x, T y)
Definition: svm.cpp:70
unsigned short PixelType
TOutputImage::PixelType OutputImagePixelType
Vector< double, TOutputImage::ImageDimension > ItkVectorType
SmartPointer< const Self > ConstPointer
#define Min(x, y)
Definition: AnnotationP.h:38
itk::ImageRegion< 3 >::SizeValueType SizeRegionType
ImageSource< TOutputImage > Superclass