18 #ifndef _FUEL_CELL_APPLICATION_CORE_BLOCK_MATRIX_APPLICATION_H_ 
   19 #define _FUEL_CELL_APPLICATION_CORE_BLOCK_MATRIX_APPLICATION_H_ 
   22 #include <deal.II/base/quadrature.h> 
   23 #include <deal.II/lac/petsc_solver.h> 
   24 #include <deal.II/lac/petsc_parallel_block_sparse_matrix.h> 
   25 #include <deal.II/lac/petsc_precondition.h> 
   26 #include <deal.II/lac/block_sparsity_pattern.h> 
   27 #include <deal.II/lac/solver_bicgstab.h> 
   28 #include <deal.II/lac/solver_cg.h> 
   29 #include <deal.II/numerics/matrix_tools.h> 
   30 #include <deal.II/base/timer.h> 
   39 #include <boost/filesystem.hpp> 
   43     namespace ApplicationCore
 
   92             boost::shared_ptr<ApplicationData>());
 
  106                                    bool triangulation_only);
 
  158             virtual void initialize (ParameterHandler& param);
 
  165             #ifdef OPENFCST_WITH_PETSC 
  226                                       const double delta = 1e-6);
 
  237             #ifdef OPENFCST_WITH_PETSC 
  238                 void residual_constraints(PETScWrappers::MPI::Vector& dst, 
const std::vector<unsigned int>& idx_list) 
const;
 
  286             #ifdef OPENFCST_WITH_PETSC 
  298           template <
typename Vector>
 
  351             #ifdef OPENFCST_WITH_PETSC 
  352                 PETScWrappers::MPI::SparseMatrix 
matrix;
 
void residual_constraints(FEVector &dst) const 
Redefinition of residual_constraints() in DoFHandler. 
 
virtual void initialize(ParameterHandler ¶m)
Initialize application. 
 
const unsigned int dim
Definition: fcst_constants.h:23
 
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:370
 
Application handling matrices and assembling linear systems of equations. 
Definition: block_matrix_application.h:74
 
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:315
 
boost::shared_ptr< Quadrature< dim-1 > > quadrature_assemble_face
Quadrature rule for matrix assembling on faces. 
Definition: block_matrix_application.h:310
 
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 "Grid Generation" in the parameter file...
 
bool print_debug
Flag specifying if the matrix and rhs should be printed at each iteration. 
Definition: block_matrix_application.h:382
 
void print_matrix_and_rhs(Vector &sys_rhs) const 
Internal routine to print matrix and rhs. 
 
BlockSparseMatrix< double > matrix
Storage for the actual matrix. 
Definition: block_matrix_application.h:354
 
bool symmetric_matrix_flag
Variable used to specify if the matrix is symmetric to the linear solver. 
Definition: block_matrix_application.h:378
 
This class is created for the objects handed to the mesh loops. 
Definition: mesh_loop_info_objects.h:544
 
bool mumps_additional_mem
Variable used for configuring MUMPS to use more memory in solve() function. 
Definition: block_matrix_application.h:374
 
void _initialize(ParameterHandler ¶m)
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:328
 
virtual void declare_parameters(ParameterHandler ¶m)
Declare parameters related to the linear system. 
 
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:321
 
bool output_system_assembling_time
If true, information about duration of the system assembling stage will be output. 
Definition: block_matrix_application.h:335
 
Timer timer
Object used to calculate the CPU and Run time for the system assembling stage. 
Definition: block_matrix_application.h:340
 
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:360
 
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:46
 
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:107
 
boost::shared_ptr< ApplicationData > data
Object for auxiliary data. 
Definition: application_base.h:348
 
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:304