DOpE
directlinearsolver.h
Go to the documentation of this file.
1 
24 #ifndef _DIRECT_LINEAR_SOLVER_H_
25 #define _DIRECT_LINEAR_SOLVER_H_
26 
27 #include <lac/vector.h>
28 #include <lac/block_sparsity_pattern.h>
29 #include <lac/block_sparse_matrix.h>
30 #include <lac/compressed_simple_sparsity_pattern.h>
31 #include <lac/solver_cg.h>
32 #include <lac/precondition.h>
33 #include <lac/full_matrix.h>
34 #include <lac/sparse_direct.h>
35 
36 #include <dofs/dof_tools.h>
37 
38 #include <numerics/vector_tools.h>
39 
40 #include <vector>
41 
42 #include "parameterreader.h"
43 
44 namespace DOpE
45 {
59  template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
61  {
62  public:
65 
66  static void declare_params(ParameterReader &param_reader);
67 
72  template<typename PROBLEM>
73  void ReInit(PROBLEM& pde);
74 
93  template<typename PROBLEM, typename INTEGRATOR>
94  void Solve(PROBLEM& pde, INTEGRATOR& integr, VECTOR &rhs, VECTOR &solution, bool force_matrix_build=false);
95 
96  protected:
97 
98  private:
99  SPARSITYPATTERN _sparsity_pattern;
100  MATRIX _matrix;
101 
102  dealii::SparseDirectUMFPACK* _A_direct;
103 
104  double _linear_global_tol, _linear_tol;
105  int _linear_maxiter;
106  };
107 
108 /*********************************Implementation************************************************/
109 
110  template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
112  {
113  }
114 
115  /******************************************************/
116 
117  template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
119  ParameterReader &/*param_reader*/)
120  {
121  _A_direct = NULL;
122  }
123 
124 /******************************************************/
125 
126 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
128 {
129  if(_A_direct != NULL)
130  {
131  delete _A_direct;
132  }
133 }
134 
135 /******************************************************/
136 
137 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
138  template<typename PROBLEM>
140 {
141  _matrix.clear();
142  pde.ComputeSparsityPattern(_sparsity_pattern);
143  _matrix.reinit(_sparsity_pattern);
144 
145  if(_A_direct != NULL)
146  {
147  delete _A_direct;
148  _A_direct= NULL;
149  }
150 }
151 
152 /******************************************************/
153 
154 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
155  template<typename PROBLEM, typename INTEGRATOR>
157  INTEGRATOR& integr,
158  VECTOR &rhs,
159  VECTOR &solution,
160  bool force_matrix_build)
161 {
162  if(force_matrix_build)
163  {
164  integr.ComputeMatrix (pde,_matrix);
165  }
166 
167  integr.ApplyNewtonBoundaryValues(pde,_matrix,rhs,solution);
168 
169  if(_A_direct == NULL)
170  {
171  _A_direct = new dealii::SparseDirectUMFPACK;
172  _A_direct->initialize(_matrix);
173  }
174  else if(force_matrix_build)
175  {
176  _A_direct->factorize(_matrix);
177  }
178 
179  dealii::Vector<double> sol;
180  sol = rhs;
181  _A_direct->solve(sol);
182  solution = sol;
183 
184  pde.GetDoFConstraints().distribute(solution);
185 
186 }
187 
188 
189 }
190 #endif
Definition: parameterreader.h:36
~DirectLinearSolverWithMatrix()
Definition: directlinearsolver.h:127
static void declare_params(ParameterReader &param_reader)
Definition: directlinearsolver.h:111
void ReInit(PROBLEM &pde)
Definition: directlinearsolver.h:139
DirectLinearSolverWithMatrix(ParameterReader &param_reader)
Definition: directlinearsolver.h:118
void Solve(PROBLEM &pde, INTEGRATOR &integr, VECTOR &rhs, VECTOR &solution, bool force_matrix_build=false)
Definition: directlinearsolver.h:156
Definition: directlinearsolver.h:60