OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
linear_solvers.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-2013 by Energy Systems Design Laboratory, University of Alberta
6 //
7 // This software is distributed under the MIT License
8 // For more information, see the README file in /doc/LICENSE
9 //
10 // - Class: linear_solvers.h
11 // - Description: This namespace contains various linear solvers and preconditioners
12 // - Developers: Valentin N. Zingan, University of Alberta
13 // Marc Secanell Gallart, University of Alberta
14 // - Id: $Id: linear_solvers.h 2605 2014-08-15 03:36:44Z secanell $
15 //
16 // ----------------------------------------------------------------------------
17 
18 #ifndef _FCST_LINEAR_SOLVERS_H_
19 #define _FCST_LINEAR_SOLVERS_H_
20 
21 #include <lac/vector.h>
22 #include <lac/block_vector.h>
23 
24 #include <lac/sparse_matrix.h>
25 #include <lac/block_sparse_matrix.h>
26 
27 #include <lac/block_matrix_array.h>
28 #include <lac/vector_memory.h>
29 
30 #include <lac/sparse_direct.h>
31 #include <lac/solver_gmres.h>
32 
33 #include <lac/precondition.h>
34 #include <lac/sparse_ilu.h>
35 
36 using namespace dealii;
37 
57 namespace LinearSolvers
58 {
59 
63  static GrowingVectorMemory< Vector<double> > vector_pool;
64  static GrowingVectorMemory< BlockVector<double> > block_vector_pool;
65 
76  {
77  public:
78 
80 
81 
85  SparseDirectUMFPACKSolver() { FcstUtilities::log << "Solving the linear system using UMFPACK" << std::endl; }
86 
91 
93 
95 
96 
106  template<typename MATRIX, typename VECTOR>
107  void solve(const MATRIX& matrix,
108  VECTOR& solution,
109  const VECTOR& right_hand_side) const
110  {
111  Vector<double> tmp;
112  tmp = right_hand_side;
113 
114  SparseDirectUMFPACK solver;
115  solver.initialize(matrix);
116  solver.solve(tmp);
117 
118  solution = tmp;
119  }
120 
122 
123  };
124 
134  {
135  public:
136 
138 
139 
143  GMRESSolver() { FcstUtilities::log << "Solving the linear system using GMRES" << std::endl; }
144 
149 
151 
153 
154 
166  template<typename MATRIX, typename VECTOR, typename PRECONDITIONER>
167  void solve(const MATRIX& matrix,
168  VECTOR& solution,
169  const VECTOR& right_hand_side,
170  const unsigned int& n_iter,
171  const double& tolerance,
172  const PRECONDITIONER& preconditioner) const
173  {
174  SolverControl solver_control(n_iter, tolerance);
175 
176  typename SolverGMRES< VECTOR >::AdditionalData additional_data(100, false, true);
177 
178  SolverGMRES< VECTOR > solver(solver_control, additional_data);
179 
180  solver.solve(matrix,
181  solution,
182  right_hand_side,
183  preconditioner);
184 
185  }
186 
198  template<typename MATRIX, typename VECTOR, typename PRECONDITIONER>
199  void solve(SolverControl& solver_control,
200  const MATRIX& matrix,
201  VECTOR& solution,
202  const VECTOR& right_hand_side,
203  const PRECONDITIONER& preconditioner) const
204  {
205  SolverGMRES< VECTOR > solver(solver_control);
206 
207  solver.solve(matrix,
208  solution,
209  right_hand_side,
210  preconditioner);
211 
212  }
214 
215  };
216 
226  {
227  public:
228 
230 
231 
235  ILUPreconditioner(const BlockSparseMatrix<double>& matrix)
236  :
237  preconditioner_matrix(matrix.n_block_rows()),
238  preconditioner_pointer(matrix.n_block_rows()),
239  preconditioner(matrix.n_block_rows())
240  {
241  FcstUtilities::log << "ILUPreconditioner is used" << std::endl;
242 
243  for(unsigned int i = 0; i < matrix.n_block_rows(); ++i)
244  {
245  preconditioner_matrix[i].initialize(matrix.block(i,i),
246  SparseILU<double>::AdditionalData());
247  preconditioner_pointer[i].set_memory(&vector_pool);
248  preconditioner_pointer[i] = &preconditioner_matrix[i];
249  }
250 
251  for(unsigned int i = 0; i < matrix.n_block_rows(); ++i)
252  for(unsigned int j = 0; j < matrix.n_block_cols(); ++j)
253  {
254  if( i == j )
255  preconditioner.enter(preconditioner_pointer[i], i, i);
256  else
257  preconditioner.enter(matrix.block(i,j), i, j);
258  }
259  }
260 
265 
267 
269 
270 
274  BlockTrianglePrecondition<double> preconditioner;
275 
276 
280  std::vector< SparseILU<double> > preconditioner_matrix;
281 
285  std::vector< PointerMatrixAux< SparseILU<double>, Vector<double> > > preconditioner_pointer;
286 
287 
289 
290  };
291 
292 } // LinearSolvers
293 
294 #endif
This class implements ILU preconditioner.
Definition: linear_solvers.h:225
SparseDirectUMFPACKSolver()
Constructor.
Definition: linear_solvers.h:85
~GMRESSolver()
Destructor.
Definition: linear_solvers.h:148
void solve(const MATRIX &matrix, VECTOR &solution, const VECTOR &right_hand_side) const
This function solves the linear system Ax = b where.
Definition: linear_solvers.h:107
This class implements GMRES solver.
Definition: linear_solvers.h:133
std::vector< PointerMatrixAux< SparseILU< double >, Vector< double > > > preconditioner_pointer
Preconditioner pointer.
Definition: linear_solvers.h:285
static GrowingVectorMemory< BlockVector< double > > block_vector_pool
Definition: linear_solvers.h:64
void solve(SolverControl &solver_control, const MATRIX &matrix, VECTOR &solution, const VECTOR &right_hand_side, const PRECONDITIONER &preconditioner) const
This routine is used when solver_control is an object already initialized.
Definition: linear_solvers.h:199
FCSTLogStream log
Object used to output data to file and, if file attached recorded to a file as well.
BlockTrianglePrecondition< double > preconditioner
Preconditioner.
Definition: linear_solvers.h:274
This class implements an interface to the sparse direct solver UMFPACK, see the link below http://www...
Definition: linear_solvers.h:75
~ILUPreconditioner()
Destructor.
Definition: linear_solvers.h:264
void solve(const MATRIX &matrix, VECTOR &solution, const VECTOR &right_hand_side, const unsigned int &n_iter, const double &tolerance, const PRECONDITIONER &preconditioner) const
This function solves the linear system Ax = b where.
Definition: linear_solvers.h:167
ILUPreconditioner(const BlockSparseMatrix< double > &matrix)
Constructor.
Definition: linear_solvers.h:235
static GrowingVectorMemory< Vector< double > > vector_pool
Memory management objects.
Definition: linear_solvers.h:63
std::vector< SparseILU< double > > preconditioner_matrix
Preconditioner matrix.
Definition: linear_solvers.h:280
GMRESSolver()
Constructor.
Definition: linear_solvers.h:143
~SparseDirectUMFPACKSolver()
Destructor.
Definition: linear_solvers.h:90