ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkFEMSolver.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 itkFEMSolver_h
20 #define itkFEMSolver_h
21 
22 #include "itkProcessObject.h"
23 #include "itkFEMObject.h"
24 
27 
28 #include "itkImage.h"
29 
30 namespace itk
31 {
32 namespace fem
33 {
68 template <unsigned int VDimension = 3>
69 class Solver : public ProcessObject
70 {
71 public:
72 
74  typedef Solver Self;
78 
80  itkNewMacro(Self);
81 
83  itkTypeMacro(Solver, ProcessObject);
84 
85  itkStaticConstMacro(FEMDimension, unsigned int, VDimension);
86  itkStaticConstMacro(MaxDimensions, unsigned int, 3);
87 
93 
101 
113 
114  /*
115  * Get/Set the Interpolation Grid Origin
116  */
117  itkSetMacro(Origin, InterpolationGridPointType);
118  itkGetMacro(Origin, InterpolationGridPointType);
119 
120  /*
121  * Get/Set the Interpolation Grid Spacing
122  */
123  itkSetMacro(Spacing, InterpolationGridSpacingType);
124  itkGetMacro(Spacing, InterpolationGridSpacingType);
125 
126  /*
127  * Get/Set the Interpolation Grid Region
128  */
129  itkSetMacro(Region, InterpolationGridRegionType);
130  itkGetMacro(Region, InterpolationGridRegionType);
131 
132  /*
133  * Get/Set the Interpolation Grid Direction
134  */
135  itkSetMacro(Direction, InterpolationGridDirectionType);
136  itkGetMacro(Direction, InterpolationGridDirectionType);
137 
139  virtual Float GetTimeStep() const;
140 
146  virtual void SetTimeStep(Float dt);
147 
149  Float GetSolution(unsigned int i, unsigned int which = 0);
150 
153  using Superclass::SetInput;
154  virtual void SetInput( FEMObjectType *fem);
155 
156  virtual void SetInput( unsigned int, FEMObjectType * fem);
157 
159 
160  FEMObjectType * GetInput(unsigned int idx);
161 
170  const Element * GetElementAtPoint(const VectorType & pt) const;
171 
173  Float GetDeformationEnergy(unsigned int SolutionIndex = 0);
174 
190 
197  {
198  return m_ls;
199  }
200 
216  const InterpolationGridPointType & bb2);
217 
220  {
223 
224  bb1.Fill(0.0);
225 
227  for( unsigned int i = 0; i < FEMDimension; i++ )
228  {
229  bb2[i] = size[i] - 1.0;
230  }
231  InitializeInterpolationGrid(size, bb1, bb2);
232  }
233 
238  const InterpolationGridPointType & origin,
239  const InterpolationGridSpacingType & spacing,
240  const InterpolationGridDirectionType & direction);
241 
253  {
255  }
256 
257 
263 
278 
279  FEMObjectType * GetOutput(unsigned int idx);
280 
281 protected:
282  Solver();
283  virtual ~Solver();
284  virtual void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
285 
288  virtual void GenerateData() ITK_OVERRIDE;
289 
290 
299  // void GenerateGFN();
300 
304  void AssembleK();
305 
312  virtual void InitializeMatrixForAssembly(unsigned int N);
313 
319  virtual void FinalizeMatrixAfterAssembly(void)
320  {
321  // Apply the boundary conditions to the K matrix
322  this->ApplyBC();
323  }
324 
332 
341 
355  void ApplyBC(int dim = 0, unsigned int matrix = 0);
356 
364  void AssembleF(int dim = 0);
365 
369  void DecomposeK();
370 
374  virtual void RunSolver();
375 
380  void UpdateDisplacements();
381 
385  void FillInterpolationGrid();
386 
391  virtual void InitializeLinearSystemWrapper();
392 
396  unsigned int m_NGFN;
397 
402  unsigned int m_NMFC;
403 
406 
411 
418 
420 
421 private:
422  Solver(const Self &); // purposely not implemented
423  void operator=(const Self &); // purposely not implemented
424 
425  /*
426  * Properties of the Interpolation Grid
427  */
432 
433 };
434 } // end namespace fem
435 } // end namespace itk
436 
437 #ifndef ITK_MANUAL_INSTANTIATION
438 #include "itkFEMSolver.hxx"
439 #endif
440 
441 #endif // #ifndef itkFEMSolver_h
void ApplyBC(int dim=0, unsigned int matrix=0)
Float GetSolution(unsigned int i, unsigned int which=0)
Superclass::RegionType RegionType
Definition: itkImage.h:140
unsigned int m_NMFC
Definition: itkFEMSolver.h:402
void InitializeInterpolationGrid(const InterpolationGridSizeType &size)
Definition: itkFEMSolver.h:219
Light weight base class for most itk classes.
InterpolationGridType::PointType InterpolationGridPointType
Definition: itkFEMSolver.h:109
Implements N-dimensional Finite element (FE) models including elements, materials, and loads.
Definition: itkFEMObject.h:75
void FillInterpolationGrid()
void SetLinearSystemWrapper(LinearSystemWrapper::Pointer ls)
virtual void FinalizeMatrixAfterAssembly(void)
Definition: itkFEMSolver.h:319
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
InterpolationGridType::SpacingType InterpolationGridSpacingType
Definition: itkFEMSolver.h:110
void InitializeInterpolationGrid(const InterpolationGridSizeType &size, const InterpolationGridPointType &bb1, const InterpolationGridPointType &bb2)
InterpolationGridType::DirectionType InterpolationGridDirectionType
Definition: itkFEMSolver.h:112
unsigned int m_NGFN
Definition: itkFEMSolver.h:396
InterpolationGridPointerType m_InterpolationGrid
Definition: itkFEMSolver.h:417
InterpolationGridSpacingType m_Spacing
Definition: itkFEMSolver.h:430
virtual void InitializeLinearSystemWrapper()
Element::VectorType VectorType
Definition: itkFEMSolver.h:96
FEMObjectType::ConstPointer FEMObjectConstPointer
Definition: itkFEMSolver.h:91
FEMObjectType::Pointer FEMObjectPointer
Definition: itkFEMSolver.h:90
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType) override
InterpolationGridDirectionType m_Direction
Definition: itkFEMSolver.h:431
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
virtual void PrintSelf(std::ostream &os, Indent indent) const override
itk::Image< Element::ConstPointer, VDimension > InterpolationGridType
Definition: itkFEMSolver.h:105
ObjectType * GetPointer() const
static const double e
The base of the natural logarithm or Euler&#39;s number
Definition: itkMath.h:45
void Fill(const ValueType &)
itk::fem::FEMObject< VDimension > FEMObjectType
Definition: itkFEMSolver.h:89
Load::ArrayType LoadArray
Definition: itkFEMSolver.h:99
SmartPointer< Self > Pointer
Definition: itkFEMSolver.h:76
InterpolationGridType::RegionType InterpolationGridRegionType
Definition: itkFEMSolver.h:108
FEMObjectPointer m_FEMObject
Definition: itkFEMSolver.h:419
static const unsigned int MaxDimensions
Definition: itkFEMSolver.h:86
virtual void AssembleElementMatrix(Element::Pointer e)
LinearSystemWrapper::Pointer m_ls
Definition: itkFEMSolver.h:405
FEMObjectType * GetOutput()
Defines all functions required by Solver class to allocate, assemble and solve a linear system of equ...
InterpolationGridType::IndexType InterpolationGridIndexType
Definition: itkFEMSolver.h:111
virtual void AssembleLandmarkContribution(Element::ConstPointer e, float)
LinearSystemWrapper::Pointer GetLinearSystemWrapper()
Definition: itkFEMSolver.h:196
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
virtual void SetTimeStep(Float dt)
Superclass::IndexType IndexType
Definition: itkImage.h:122
vnl_vector< Float > VectorType
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
InterpolationGridRegionType m_Region
Definition: itkFEMSolver.h:428
Element::ArrayType ElementArray
Definition: itkFEMSolver.h:98
InterpolationGridType::Pointer InterpolationGridPointerType
Definition: itkFEMSolver.h:106
InterpolationGridPointType m_Origin
Definition: itkFEMSolver.h:429
FEM solver used to generate a solution for a FE formulation.
Definition: itkFEMSolver.h:69
LinearSystemWrapperVNL m_lsVNL
Definition: itkFEMSolver.h:410
A region represents some portion or piece of data.
Definition: itkRegion.h:64
virtual void InitializeMatrixForAssembly(unsigned int N)
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx)
Make a DataObject of the correct type to used as the specified output.
Element::Node::ArrayType NodeArray
Definition: itkFEMSolver.h:97
virtual void SetInput(FEMObjectType *fem)
void operator=(const Self &)
Array for FEMP objects.
Definition: itkFEMPArray.h:40
ProcessObject Superclass
Definition: itkFEMSolver.h:75
void AssembleF(int dim=0)
Abstract base element class.
InterpolationGridType::SizeType InterpolationGridSizeType
Definition: itkFEMSolver.h:107
static const unsigned int FEMDimension
Definition: itkFEMSolver.h:85
void UpdateDisplacements()
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Definition: itkFEMSolver.h:260
DataObject::Pointer DataObjectPointer
Definition: itkFEMSolver.h:92
virtual Float GetTimeStep() const
const InterpolationGridType * GetInterpolationGrid(void) const
Definition: itkFEMSolver.h:252
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Float GetDeformationEnergy(unsigned int SolutionIndex=0)
FEMObjectType * GetInput()
virtual void GenerateData() override
const Element * GetElementAtPoint(const VectorType &pt) const
Element::Float Float
Definition: itkFEMSolver.h:95
virtual void RunSolver()
SmartPointer< const Self > ConstPointer
Definition: itkFEMSolver.h:77
LinearSystemWrapper class that uses VNL numeric library functions to define a sparse linear system of...
Templated n-dimensional image class.
Definition: itkImage.h:75
Material::ArrayType MaterialArray
Definition: itkFEMSolver.h:100