OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
block_matrix_application.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-2009 by Guido Kanschat
6 // Copyright (C) 2006-2014 by Energy Systems Design Laboratory, University of Alberta
7 //
8 // This software is distributed under the MIT License
9 // For more information, see the README file in /doc/LICENSE
10 //
11 // - Class: block_matrix_application.h
12 // - Description:
13 // - Developers: Guido Kanschat, Texas A&M University
14 // Marc Secanell, University of Alberta
15 // - Id: $Id: block_matrix_application.h 2605 2014-08-15 03:36:44Z secanell $
16 //
17 // ----------------------------------------------------------------------------
18 
19 #ifndef _FUEL_CELL_APPLICATION_CORE_BLOCK_MATRIX_APPLICATION_H_
20 #define _FUEL_CELL_APPLICATION_CORE_BLOCK_MATRIX_APPLICATION_H_
21 
22 #include <base/quadrature.h>
23 #include <lac/full_matrix.h>
24 #include <lac/block_sparse_matrix.h>
25 #include <lac/block_sparsity_pattern.h>
26 #include <lac/solver_bicgstab.h>
27 #include <lac/solver_cg.h>
28 #include <lac/vector.h>
29 #include <grid/tria.h>
30 #include <grid/tria_iterator.h>
31 #include <grid/tria_accessor.h>
32 #include <dofs/dof_tools.h>
33 #include <numerics/matrix_tools.h>
34 #include <base/parameter_handler.h>
35 #include <base/quadrature_lib.h>
36 #include <dofs/dof_handler.h>
37 #include <fe/fe.h>
38 #include <fe/fe_values.h>
39 
40 
41 #include <lac/petsc_solver.h>
42 #include <lac/petsc_parallel_block_sparse_matrix.h>
43 #include <lac/petsc_precondition.h>
44 
45 #include "system_management.h"
46 #include "dof_application.h"
47 #include "matrix_block.h"
48 #include "linear_solvers.h"
49 
50 namespace FuelCell
51 {
52  namespace ApplicationCore
53  {
82  template <int dim>
84  public DoFApplication<dim>
85  {
86  public:
88 
89 
100  BlockMatrixApplication(boost::shared_ptr<ApplicationData> data =
101  boost::shared_ptr<ApplicationData>());
102 
115  bool triangulation_only);
116 
124  void _initialize(ParameterHandler& param);
125 
161  virtual void declare_parameters(ParameterHandler& param);
165  virtual void initialize (ParameterHandler& param);
167 
172  #ifdef OPENFCST_WITH_PETSC
173  void PETSc_assemble(const FEVectors&);
174  #else
175  void serial_assemble(const FEVectors&);
176  #endif
177 
178 
179 
180 
186  virtual void post_cell_assemble();
187 
189 
190 
195  void remesh_matrices();
200  virtual void remesh();
201 
208  void assemble(const FEVectors&);
209 
210 
232  void assemble_numerically(const FEVectors& src,
233  const double delta = 1e-6);
234 
242  void residual_constraints(FEVector& dst) const;
243 
244  #ifdef OPENFCST_WITH_PETSC
245  void residual_constraints(PETScWrappers::MPI::Vector& dst, const std::vector<unsigned int>& idx_list) const;
246  #endif
247 
249 
251 
257  virtual void cell_matrix(MatrixVector& cell_matrices,
258  const typename DoFApplication<dim>::CellInfo& cell);
259 
264  virtual void bdry_matrix(MatrixVector& face_matrices,
265  const typename DoFApplication<dim>::FaceInfo& face);
266 
271  virtual void face_matrix(MatrixVector& matrices11,
272  MatrixVector& matrices12,
273  MatrixVector& matrices21,
274  MatrixVector& matrices22,
275  const typename DoFApplication<dim>::FaceInfo& face1,
276  const typename DoFApplication<dim>::FaceInfo& face2);
277 
285 
287  virtual void dirichlet_bc(std::map<unsigned int, double>& boundary_values) const;
289 
290  virtual void solve(FEVector& dst,
291  const FEVectors& src);
292 
293  #ifdef OPENFCST_WITH_PETSC
294  void PETSc_solve(FuelCell::ApplicationCore::FEVector system_rhs,
295  FEVector& solution, const FEVectors& src);
296  #else
298  FEVector& solution);
299  #endif
300 
301 
302  protected:
307  boost::shared_ptr<Quadrature<dim> > quadrature_assemble_cell;
308 
313  boost::shared_ptr<Quadrature<dim-1> > quadrature_assemble_face;
314 
319 
320 
324  std::map<unsigned int, double> boundary_values;
325 
326  private:
327 
331  BlockSparsityPattern sparsities;
332 
333 
334  protected:
336 
337 
342  #ifdef OPENFCST_WITH_PETSC
343  PETScWrappers::MPI::SparseMatrix matrix;
344  #else
345  BlockSparseMatrix<double> matrix;
346  #endif
347 
351  SolverControl solver_control;
352 
353 
355 
357 
358 
361  std::string solver_name;
371  };
372 
373 }
374 
375 }
376 
377 #endif
void residual_constraints(FEVector &dst) const
Redefinition of residual_constraints() in DoFHandler.
virtual void initialize(ParameterHandler &param)
Initialize application.
const unsigned int dim
Definition: fcst_constants.h:24
void assemble(const FEVectors &)
Loop over all cells and assemble the system #matrices by using the local matrix integration functions...
bool assemble_numerically_flag
Variable used to select if assembly should be done analytically or numerically.
Definition: block_matrix_application.h:365
Application handling matrices and assembling linear systems of equations.
Definition: block_matrix_application.h:83
virtual void dirichlet_bc(std::map< unsigned int, double > &boundary_values) const
bool repair_diagonal
Bool determining whether or not to repair diagonal before solving.
Definition: block_matrix_application.h:318
boost::shared_ptr< Quadrature< dim-1 > > quadrature_assemble_face
Quadrature rule for matrix assembling on faces.
Definition: block_matrix_application.h:313
BlockMatrixApplication(boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >())
Constructor for an object, owning its own mesh and dof handler.
void serial_assemble(const FEVectors &)
serial and PETSc assemble functions called by assemble().
virtual void post_cell_assemble()
A call back function for assemble(), called after all cell matrices have been entered, but before the face matrices are computed.
void assemble_numerically(const FEVectors &src, const double delta=1e-6)
Compute the Jacobian of the system of equations, by using forward differences.
virtual void remesh()
Refine grid accordingly to the Refinement options under &quot;Grid Generation&quot; in the parameter file...
BlockSparseMatrix< double > matrix
Storage for the actual matrix.
Definition: block_matrix_application.h:345
This class is created for the objects handed to the mesh loops.
Definition: mesh_loop_info_objects.h:625
bool mumps_additional_mem
Variable used for configuring MUMPS to use more memory in solve() function.
Definition: block_matrix_application.h:369
std::string solver_name
Variable used to select the appropriate linear solver:
Definition: block_matrix_application.h:361
void _initialize(ParameterHandler &param)
Initialize data of this class.
std::vector< MatrixBlock< FullMatrix< double > > > MatrixVector
The matrix vector used in the mesh loops.
Definition: matrix_block.h:102
BlockSparsityPattern sparsities
Sparsity patterns.
Definition: block_matrix_application.h:331
virtual void declare_parameters(ParameterHandler &param)
Declare parameters related to matrices.
std::map< unsigned int, double > boundary_values
Variable to store boundary values, so they only need to be computed once per mesh refinement...
Definition: block_matrix_application.h:324
virtual void solve(FEVector &dst, const FEVectors &src)
Solve the system assembled with right hand side in FEVectors src and return the result in FEVector ds...
void remesh_matrices()
Initialize sparsity patterns and matrices for the new mesh.
SolverControl solver_control
Solver control object.
Definition: block_matrix_application.h:351
virtual void face_matrix(MatrixVector &matrices11, MatrixVector &matrices12, MatrixVector &matrices21, MatrixVector &matrices22, const typename DoFApplication< dim >::FaceInfo &face1, const typename DoFApplication< dim >::FaceInfo &face2)
Integration of local bilinear form.
virtual void bdry_matrix(MatrixVector &face_matrices, const typename DoFApplication< dim >::FaceInfo &face)
Integration of local bilinear form.
void serial_solve(FuelCell::ApplicationCore::FEVector system_rhs, FEVector &solution)
BlockVector< double > FEVector
The vector class used by applications.
Definition: application_data.h:39
The data type used in function calls of Application.
Definition: fe_vectors.h:59
Base class for all linear applications, i.e., all applications requiring Triangulation and DoFHandler...
Definition: dof_application.h:122
boost::shared_ptr< ApplicationData > data
Object for auxiliary data.
Definition: application_base.h:342
virtual void cell_matrix(MatrixVector &cell_matrices, const typename DoFApplication< dim >::CellInfo &cell)
Integration of local bilinear form.
boost::shared_ptr< Quadrature< dim > > quadrature_assemble_cell
Quadrature rule for matrix assembling on cells.
Definition: block_matrix_application.h:307