Statismo  0.10.1
 All Classes Namespaces Functions Typedefs
itkStandardImageRepresenter.h
1 /*
2  * This file is part of the statismo library.
3  *
4  * Author: Marcel Luethi (marcel.luethi@unibas.ch)
5  *
6  * Copyright (c) 2011 University of Basel
7  * All rights reserved.
8  *
9  * Redistribution and use in source and binary forms, with or without
10  * modification, are permitted provided that the following conditions
11  * are met:
12  *
13  * Redistributions of source code must retain the above copyright notice,
14  * this list of conditions and the following disclaimer.
15  *
16  * Redistributions in binary form must reproduce the above copyright
17  * notice, this list of conditions and the following disclaimer in the
18  * documentation and/or other materials provided with the distribution.
19  *
20  * Neither the name of the project's author nor the names of its
21  * contributors may be used to endorse or promote products derived from
22  * this software without specific prior written permission.
23  *
24  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
28  * HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
29  * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
30  * TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
31  * PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
32  * LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
33  * NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
34  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
35  *
36  */
37 
38 #ifndef ITK_STANDARDIMAGE_REPRESENTER_H_
39 #define ITK_STANDARDIMAGE_REPRESENTER_H_
40 
41 #include "statismoITKConfig.h" // this needs to be the first include
42 
43 #include <H5Cpp.h>
44 
45 #include <itkObject.h>
46 #include <itkImage.h>
47 
48 #include "itkPixelConversionTraits.h"
49 #include "itkStandardImageRepresenterTraits.h"
50 
51 #include "CommonTypes.h"
52 #include "Representer.h"
53 
54 namespace itk {
55 
62 template<class TPixel, unsigned ImageDimension>
63 class StandardImageRepresenter: public Object, public statismo::Representer<
64  itk::Image<TPixel, ImageDimension> > {
65  public:
66 
67  /* Standard class typedefs. */
69  typedef Object Superclass;
70  typedef SmartPointer<Self> Pointer;
71  typedef SmartPointer<const Self> ConstPointer;
72 
75 
78 
79  typedef itk::Image<TPixel, ImageDimension> ImageType;
80 
81  typedef typename statismo::Representer<ImageType> RepresenterBaseType;
83  typedef typename RepresenterBaseType::PointType PointType;
84  typedef typename RepresenterBaseType::ValueType ValueType;
85  typedef typename RepresenterBaseType::DatasetType DatasetType;
86  typedef typename RepresenterBaseType::DatasetPointerType DatasetPointerType;
87  typedef typename RepresenterBaseType::DatasetConstPointerType DatasetConstPointerType;
88 
89  static StandardImageRepresenter* Create() {
90  return new StandardImageRepresenter();
91  }
92  void Load(const H5::Group& fg);
93  StandardImageRepresenter* Clone() const;
94 
96  virtual ~StandardImageRepresenter();
97 
98  unsigned GetDimensions() const {
99  return PixelConversionTrait<TPixel>::GetPixelDimension();
100  }
101  std::string GetName() const {
102  return "itkStandardImageRepresenter";
103  }
104  typename RepresenterBaseType::RepresenterDataType GetType() const {
105  return RepresenterBaseType::IMAGE;
106  }
107 
108  const DomainType& GetDomain() const {
109  return m_domain;
110  }
111  std::string GetVersion() const {
112  return "0.1";
113  }
114 
116  DatasetConstPointerType GetReference() const {
117  return m_reference;
118  }
119 
120 
122  void SetReference(ImageType* ds);
123 
128  statismo::VectorType PointToVector(const PointType& pt) const;
129  DatasetPointerType DatasetToSample(DatasetConstPointerType ds) const;
130  statismo::VectorType SampleToSampleVector(DatasetConstPointerType sample) const;
131  DatasetPointerType SampleVectorToSample(
132  const statismo::VectorType& sample) const;
133 
134  ValueType PointSampleFromSample(DatasetConstPointerType sample,
135  unsigned ptid) const;
136  ValueType PointSampleVectorToPointSample(
137  const statismo::VectorType& pointSample) const;
138  statismo::VectorType PointSampleToPointSampleVector(
139  const ValueType& v) const;
140 
141  void Save(const H5::Group& fg) const;
142  virtual unsigned GetPointIdForPoint(const PointType& point) const;
143 
144  unsigned GetNumberOfPoints() const;
145 
146  void Delete() const {
147  this->UnRegister();
148  }
149 
150 
151  void DeleteDataset(DatasetConstPointerType d) const {}
152  DatasetPointerType CloneDataset(DatasetConstPointerType d) const;
153 
154  private:
155 
156  typename ImageType::Pointer LoadRef(const H5::Group& fg) const;
157  typename ImageType::Pointer LoadRefLegacy(const H5::Group& fg) const;
158 
159  DatasetConstPointerType m_reference;
160  DomainType m_domain;
161 };
162 
163 } // namespace itk
164 
165 #include "itkStandardImageRepresenter.hxx"
166 
167 #endif /* itkStandardImageREPRESENTER_H_ */
A representer for scalar and vector valued images.
Definition: itkStandardImageRepresenter.h:63
statismo::VectorType PointToVector(const PointType &pt) const
Definition: itkStandardImageRepresenter.hxx:221
itkTypeMacro(StandardImageRepresenter, Object)
void SetReference(ImageType *ds)
Definition: itkStandardImageRepresenter.hxx:202
DatasetConstPointerType GetReference() const
return the reference used in the representer
Definition: itkStandardImageRepresenter.h:116
Definition: Domain.h:51