ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkFEMSolverCrankNicolson.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 itkFEMSolverCrankNicolson_h
20 #define itkFEMSolverCrankNicolson_h
21 
22 #include "itkFEMSolver.h"
23 #include "itkFEMElementBase.h"
24 #include "itkFEMMaterialBase.h"
25 #include "itkFEMLoadBase.h"
27 
28 #include "vnl/vnl_sparse_matrix.h"
29 #include "vnl/vnl_matrix.h"
30 #include "vnl/vnl_vector.h"
31 #include "vnl/algo/vnl_svd.h"
32 #include "vnl/algo/vnl_cholesky.h"
33 #include <vnl/vnl_sparse_matrix_linear_system.h>
34 #include <cmath>
35 
36 namespace itk
37 {
38 namespace fem
39 {
70 
71 template <unsigned int TDimension = 3>
72 class ITK_TEMPLATE_EXPORT SolverCrankNicolson : public Solver<TDimension>
73 {
74 public:
79 
81  itkNewMacro(Self);
82 
85 
87 
89  itkSetMacro(UseMassMatrix, bool);
90  itkGetMacro(UseMassMatrix, bool);
92 
94  itkGetConstMacro(Iterations, unsigned int);
95 
102  void ResetIterations(void)
103  {
104  m_Iterations = 0;
105  }
106 
111  void AddToDisplacements(Float optimum = 1.0);
112 
113  void AverageLastTwoDisplacements(Float t = 0.5);
114 
115  void ZeroVector(int which = 0);
116 
117  void PrintDisplacements();
118 
119  void PrintForce();
120 
122  itkGetMacro(TotalSolutionIndex, unsigned int);
123 
125  itkGetMacro(SolutionTMinus1Index, unsigned int);
126 
128  itkSetMacro(Alpha, Float);
129  itkGetMacro(Alpha, Float);
131 
133  itkSetMacro(Rho, Float);
134  itkGetMacro(Rho, Float);
136 
138  virtual Float GetTimeStep(void) const ITK_OVERRIDE
139  {
140  return m_TimeStep;
141  }
142 
148  virtual void SetTimeStep(Float dt) ITK_OVERRIDE
149  {
150  m_TimeStep = dt;
151  }
152 
156  void RecomputeForceVector(unsigned int index);
157 
160  void FindBracketingTriplet(Float *a, Float *b, Float *c);
161 
165  Float GoldenSection(Float tol = 0.01, unsigned int MaxIters = 25);
166 
167  /* Brents method from Numerical Recipes. */
168  Float BrentsMethod(Float tol = 0.01, unsigned int MaxIters = 25);
169 
170  Float EvaluateResidual(Float t = 1.0);
171 
172  Float GetDeformationEnergy(Float t = 1.0);
173 
174  inline Float GSSign(Float a, Float b)
175  {
176  return b > 0.0 ? std::fabs(a) : -1. * std::fabs(a);
177  }
178  inline Float GSMax(Float a, Float b)
179  {
180  return a > b ? a : b;
181  }
182 
183  void SetEnergyToMin(Float xmin);
184 
186  {
187  return this->m_LinearSystem;
188  }
189 
191  {
192  return m_CurrentMaxSolution;
193  }
194 
197  void PrintMinMaxOfSolution();
198 
199 protected:
200 
206  ~SolverCrankNicolson() ITK_OVERRIDE {}
207 
210  void GenerateData() ITK_OVERRIDE;
211 
215  virtual void RunSolver(void) ITK_OVERRIDE;
216 
220  void InitializeForSolution();
221 
227  void AssembleKandM();
228 
236  void AssembleFforTimeStep(int dim = 0);
237 
238  Float m_TimeStep;
239  Float m_Rho;
240  Float m_Alpha;
241  Float m_CurrentMaxSolution;
242 
243  bool m_UseMassMatrix;
244  unsigned int m_Iterations;
245 
246  unsigned int m_ForceTIndex;
247  unsigned int m_ForceTotalIndex;
248  unsigned int m_ForceTMinus1Index;
249  unsigned int m_SolutionTIndex;
250  unsigned int m_SolutionTMinus1Index;
251  unsigned int m_SolutionVectorTMinus1Index;
252  unsigned int m_TotalSolutionIndex;
253  unsigned int m_DifferenceMatrixIndex;
254  unsigned int m_SumMatrixIndex;
255  unsigned int m_DiffMatrixBySolutionTMinus1Index;
256 
257 private:
258  ITK_DISALLOW_COPY_AND_ASSIGN(SolverCrankNicolson);
259 
260 };
261 }
262 } // end namespace itk::fem
263 
264 #ifndef ITK_MANUAL_INSTANTIATION
265 #include "itkFEMSolverCrankNicolson.hxx"
266 #endif
267 
268 #endif // #ifndef itkFEMSolverCrankNicolson_h
virtual void SetTimeStep(Float dt) ITK_OVERRIDE
Light weight base class for most itk classes.
FEM Solver for time dependent problems; uses Crank-Nicolson implicit discretization scheme...
SmartPointer< const Self > ConstPointer
Defines all functions required by Solver class to allocate, assemble and solve a linear system of equ...
virtual Float GetTimeStep(void) const ITK_OVERRIDE
FEM solver used to generate a solution for a FE formulation.
Definition: itkFEMSolver.h:73