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  };
105 
106 /*********************************Implementation************************************************/
107 
108  template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
110  {
111  }
112 
113  /******************************************************/
114 
115  template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
117  ParameterReader &/*param_reader*/)
118  {
119  A_direct_ = NULL;
120  }
121 
122 /******************************************************/
123 
124 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
126 {
127  if(A_direct_ != NULL)
128  {
129  delete A_direct_;
130  }
131 }
132 
133 /******************************************************/
134 
135 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
136  template<typename PROBLEM>
138 {
139  matrix_.clear();
140  pde.ComputeSparsityPattern(sparsity_pattern_);
141  matrix_.reinit(sparsity_pattern_);
142 
143  if(A_direct_ != NULL)
144  {
145  delete A_direct_;
146  A_direct_= NULL;
147  }
148 }
149 
150 /******************************************************/
151 
152 template <typename SPARSITYPATTERN, typename MATRIX, typename VECTOR>
153  template<typename PROBLEM, typename INTEGRATOR>
155  INTEGRATOR& integr,
156  VECTOR &rhs,
157  VECTOR &solution,
158  bool force_matrix_build)
159 {
160  if(force_matrix_build)
161  {
162  integr.ComputeMatrix (pde,matrix_);
163  }
164 
165  integr.ApplyNewtonBoundaryValues(pde,matrix_,rhs,solution);
166 
167  if(A_direct_ == NULL)
168  {
169  A_direct_ = new dealii::SparseDirectUMFPACK;
170  A_direct_->initialize(matrix_);
171  }
172  else if(force_matrix_build)
173  {
174  A_direct_->factorize(matrix_);
175  }
176 
177  dealii::Vector<double> sol;
178  sol = rhs;
179  A_direct_->solve(sol);
180  solution = sol;
181 
182  pde.GetDoFConstraints().distribute(solution);
183 
184 }
185 
186 
187 }
188 #endif
Definition: parameterreader.h:36
~DirectLinearSolverWithMatrix()
Definition: directlinearsolver.h:125
static void declare_params(ParameterReader &param_reader)
Definition: directlinearsolver.h:109
void ReInit(PROBLEM &pde)
Definition: directlinearsolver.h:137
DirectLinearSolverWithMatrix(ParameterReader &param_reader)
Definition: directlinearsolver.h:116
void Solve(PROBLEM &pde, INTEGRATOR &integr, VECTOR &rhs, VECTOR &solution, bool force_matrix_build=false)
Definition: directlinearsolver.h:154
Definition: directlinearsolver.h:60