ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkFEMScatteredDataPointSetToImageFilter.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 #ifndef __itkFEMScatteredDataPointSetToImageFilter_h
20 #define __itkFEMScatteredDataPointSetToImageFilter_h
21 
23 #include "itkVectorContainer.h"
24 #include "vnl/vnl_matrix.h"
25 
26 #include "itkMesh.h"
27 #include "itkFEMObject.h"
33 #include "itkFEMRobustSolver.h"
36 #include "itkTriangleCell.h"
37 #include "itkTetrahedronCell.h"
38 #include "itkQuadrilateralCell.h"
39 #include "itkHexahedronCell.h"
40 
41 namespace itk
42 {
43 namespace fem
44 {
126 template<typename TInputPointSet, typename TInputMesh, typename TOutputImage, typename TInputConfidencePointSet, typename TInputTensorPointSet>
128  public PointSetToImageFilter< TInputPointSet, TOutputImage >
129 {
130 public:
135 
137  itkNewMacro( Self );
138 
140  itkStaticConstMacro( ImageDimension, unsigned int, TOutputImage::ImageDimension );
141 
143  typedef TInputPointSet PointSetType;
144  typedef typename PointSetType::PointType PointType;
145  typedef typename PointSetType::PointsContainer PointsContainer;
146  typedef typename PointsContainer::ConstIterator PointsIterator;
147  typedef typename PointSetType::PixelType PointDataType;
148  typedef typename PointSetType::PointDataContainer PointDataContainerType;
149  typedef typename PointDataContainerType::ConstIterator PointDataIterator;
150 
152  typedef TInputConfidencePointSet ConfidencePointSetType;
153  typedef typename ConfidencePointSetType::PointsContainer::ConstIterator ConfidencePointsIterator;
154  typedef typename ConfidencePointSetType::PixelType ConfidencePointDataType;
155  typedef typename ConfidencePointSetType::PointDataContainer ConfidencePointDataContainerType;
156 
158  typedef TInputTensorPointSet TensorPointSetType;
159  typedef typename TensorPointSetType::PointsContainer::ConstIterator TensorPointsIterator;
160  typedef typename TensorPointSetType::PixelType TensorPointDataType;
161  typedef typename TensorPointSetType::PointDataContainer TensorPointDataContainerType;
162  typedef typename TensorPointDataContainerType::Iterator TensorPointDataIterator;
163 
165  typedef TInputMesh MeshType;
166  typedef typename MeshType::CellType CellType;
167  typedef typename CellType::CellAutoPointer CellAutoPointer;
168  typedef typename MeshType::CellsContainer CellsContainer;
169  typedef typename CellsContainer::ConstIterator CellIterator;
170 
175  typedef typename CellType::PointIdIterator PointIdIterator;
176 
178  typedef TOutputImage ImageType;
179  typedef typename ImageType::PixelType PixelType;
180  typedef typename ImageType::RegionType RegionType;
181  typedef typename ImageType::SizeType SizeType;
182  typedef typename ImageType::IndexType IndexType;
183  typedef typename ImageType::SpacingType SpacingType;
185 
187 
190 
193 
199 
202 
205 
209 
213 
219 
221  itkSetConstObjectMacro(ConfidencePointSet, ConfidencePointSetType);
222 
223  itkSetConstObjectMacro(TensorPointSet, TensorPointSetType);
224 
225  itkSetObjectMacro(Mesh, MeshType);
226  itkGetModifiableObjectMacro(Mesh, MeshType);
227 
228  itkSetObjectMacro(FEMSolver, FEMSolverType);
229  itkGetModifiableObjectMacro(FEMSolver, FEMSolverType);
230 
232  itkGetConstReferenceMacro(PixelsPerElement, ContinuousIndexType);
233  itkSetMacro(PixelsPerElement, ContinuousIndexType);
235 
237  void SetElementSpacing(const SpacingType & elementSpacing);
238  itkGetConstReferenceMacro(SpacingPerElement, SpacingType);
240 
242  itkGetConstReferenceMacro(NumberOfElements, SizeType);
243 
244 protected:
245 
248 
251 
254 
257 
259  void InitializeFEMObject(FEMObjectType * femObject);
260 
262  void InitializeMaterials(FEMObjectType * femObject);
263 
265  void InitializeNodes(FEMObjectType * femObject);
266 
268  void InitializeElements(FEMObjectType * femObject);
269 
271  void InitializeLoads(FEMObjectType * femObject);
272 
274  void GenerateData();
275 
277 
278  void PrintSelf(std::ostream & os, Indent indent) const;
279 
280 private:
281 
282  FEMScatteredDataPointSetToImageFilter( const Self & ); // purposely not implemented.
283  void operator=( const Self & ); // purposely not implemented
284 
288  typename MeshType::Pointer m_Mesh;
289 
290  typename ConfidencePointSetType::ConstPointer m_ConfidencePointSet;
291  typename TensorPointSetType::ConstPointer m_TensorPointSet;
292 
297 
300 };
301 
302 }// end namespace fem
303 }// end namespace itk
304 
305 #ifndef ITK_MANUAL_INSTANTIATION
306 #include "itkFEMScatteredDataPointSetToImageFilter.hxx"
307 #endif
308 
309 #endif
TensorPointSetType::PointsContainer::ConstIterator TensorPointsIterator
Implements N-dimensional Finite element (FE) models including elements, materials, and loads.
Definition: itkFEMObject.h:74
void InitializeElements(FEMObjectType *femObject)
void InitializeNodes(FEMObjectType *femObject)
ConfidencePointSetType::PointDataContainer ConfidencePointDataContainerType
void SetElementSpacing(const SpacingType &elementSpacing)
void PrintSelf(std::ostream &os, Indent indent) const
Generate a rectilinar mesh from an image. The result is stored in a FEMObject.
Class that stores information required to define a node.
Implements the N-dimensional mesh structure.
Definition: itkMesh.h:108
Base class for all process objects that output image data.
ImageToRectilinearFEMObjectFilter< ImageType > ImageToRectilinearFEMObjectFilterType
PointSetToImageFilter< TInputPointSet, TOutputImage > Superclass
vnl_matrix< Float > MatrixType
Scattered data approximation to interpolation in the presence of outliers.
void InitializeMaterials(FEMObjectType *femObject)
Linear elasticity material class.
TOutputImage::SpacingType SpacingType
3-noded finite element class in 2D space for linear elasticity problem.
ConfidencePointSetType::PointsContainer::ConstIterator ConfidencePointsIterator
Represents a hexahedron for a Mesh.
ContinuousIndex< SpacePrecisionType, ImageDimension > ContinuousIndexType
vnl_vector< Float > VectorType
void InitializeLoads(FEMObjectType *femObject)
4-noded finite element class in 2D space for linear elasticity problemThe ordering of the nodes is co...
TetrahedronCell represents a tetrahedron for a Mesh.
void InitializeFEMObject(FEMObjectType *femObject)
4-noded finite element class in 3D space for linear elasticity problem
Base class for filters that take a PointSet as input and produce an image as output. By default, if the user does not specify the size of the output image, the maximum size of the point-set&#39;s bounding box is used.
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Represents a quadrilateral for a Mesh.
8-noded finite element class in 3D space. The constitutive equation used is from linear elasticity th...