OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Types | Public Member Functions | Protected Attributes | Private Attributes | List of all members
AppFrame::MatrixApplication< dim > Class Template Reference

Application handling matrices and assembling linear systems of equations. More...

#include <matrix_application.h>

Inheritance diagram for AppFrame::MatrixApplication< dim >:
Inheritance graph
[legend]
Collaboration diagram for AppFrame::MatrixApplication< dim >:
Collaboration graph
[legend]

Classes

class  LocalMatrixIntegrator
 

Public Types

typedef MatrixBlock
< SparseMatrix< double > > 
MatrixType
 The type of matrices stored.
 
typedef MatrixBlock
< MGLevelObject< SparseMatrix
< double > > > 
MGMatrixType
 The type of multigrid matrices stored.
 
typedef
MeshWorker::Assembler::MatrixLocalBlocksToGlobalBlocks
< SparseMatrix< double > > 
Assembler
 The assembler object.
 
typedef
MeshWorker::Assembler::MGMatrixLocalBlocksToGlobalBlocks
< SparseMatrix< double > > 
MGAssembler
 The multigrid assembler object.
 
- Public Types inherited from AppFrame::DoFApplication< dim >
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEValuesBase< dim > > 
CellInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 
typedef
MeshWorker::InfoObjects::IntegrationInfo
< dim, FEFaceValuesBase< dim > > 
FaceInfo
 This is a redefinition of CellInfo in order to call MeshWorker objects instead of the previous implementation in AppFrame.
 

Public Member Functions

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.
 
Constructors and Initialization
 MatrixApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >(), bool multigrid=false)
 Constructor for an object, owning its own mesh and dof handler.
 
 MatrixApplication (bool multigrid)
 Constructor for an object, owning its own mesh and dof handler and having its own data.
 
 MatrixApplication (DoFApplication< dim > &, bool triangulation_only, bool multigrid=false)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
void _initialize (ParameterHandler &param)
 Initialize data of this class.
 
virtual void declare_parameters (ParameterHandler &param)
 Declare parameters related to matrices.
 
virtual void initialize (ParameterHandler &param)
 The virtual initialize function to be replaced by derived applications.
 
Data management functions
void add_matrix (unsigned int row, unsigned int col)
 Add a new matrix to the list of matrices.
 
void add_mg_matrix (unsigned int row, unsigned int col)
 Add a new multigrid matrix to mg_matrices.
 
void add_vector_for_cell_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the cell_matrix() function via CellInfo::values and CellInfo::derivatives.
 
void add_vector_for_flux_matrix (std::string name, unsigned int block, unsigned int components, bool values, bool derivatives)
 Add a named vector from data to be handed to the bdry_matrix() and face_matrix() functions via FaceInfo::values and FaceInfo::derivatives.
 
ApplicationBase interface
void remesh_matrices ()
 Initialize sparsity patterns and matrices for the new mesh.
 
virtual void remesh ()
 Refine grid accordingly to the Refinement options under "Grid Generation" in the parameter file.
 
void assemble (const FEVectors &)
 Loop over all cells and assemble the system matrices by using the local matrix integration functions cell_matrix(), bdry_matrix() and face_matrix() provided by the derived class.
 
Local integrators
virtual void cell_matrix (MatrixVector &cell_matrices, const typename DoFApplication< dim >::CellInfo &cell)
 Integration of local bilinear form.
 
virtual void bdry_matrix (MatrixVector &face_matrices, const typename DoFApplication< dim >::FaceInfo &face)
 Integration of local bilinear form.
 
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 mg_cell_matrix (MatrixVector &cell_matrices, const typename DoFApplication< dim >::CellInfo &cell)
 Integration of local bilinear form.
 
virtual void mg_bdry_matrix (MatrixVector &face_matrices, const typename DoFApplication< dim >::FaceInfo &face)
 Integration of local bilinear form.
 
virtual void mg_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.
 
- Public Member Functions inherited from AppFrame::DoFApplication< dim >
 DoFApplication (boost::shared_ptr< ApplicationData > data=boost::shared_ptr< ApplicationData >(), bool multigrid=false)
 Constructor for an object, owning its own mesh and dof handler.
 
 DoFApplication (bool multigrid)
 Constructor for an object owning its own mesh and dof handler and creating new ApplicationData.
 
 DoFApplication (DoFApplication< dim > &, bool triangulation_only, bool multigrid=false)
 Constructor for an object, borrowing mesh and dof handler from another object.
 
 ~DoFApplication ()
 Destructor which deletes owned objects.
 
void _initialize (ParameterHandler &param)
 Initialize from parameter values.
 
void add_vector_for_transfer (FEVector *)
 Add the vector to be transfered from one mesh to the next.
 
void delete_vector_for_transfer ()
 Delete the vector to be transfered from one mesh to the next.
 
virtual void remesh_dofs ()
 Initialize dof handler, count the dofs in each block and renumber the dofs.
 
virtual void init_vector (FEVector &dst) const
 Reinitialize the BlockVector dst such that it contains block_size.size() blocks.
 
virtual double evaluate (const FEVectors &)
 Evaluate a functional during postprocessing such as drag or lift on an aerodynamics problem.
 
virtual double estimate (const FEVectors &src)
 Estimate the error.
 
virtual double residual (FEVector &dst, const FEVectors &src, bool apply_boundaries=true)
 Loop over all cells to compute the residual.
 
virtual void residual_constraints (FEVector &dst) const
 Apply boundary conditions and hanging node constraints to a residual vector after it has been computed by residual().
 
unsigned int memory_consumption () const
 Compute the amount of memory needed by this object.
 
virtual double cell_estimate (const CellInfo &)
 Integration of estimate inside a cell.
 
virtual double bdry_estimate (const FaceInfo &)
 Integration of estimate on a boundary face.
 
virtual double face_estimate (const FaceInfo &, const FaceInfo &)
 Integration of estimate on an interior face.
 
virtual void cell_residual (FEVector &cell_vector, const CellInfo &cell)
 Integration of local linear form.
 
virtual void bdry_residual (FEVector &face_vector, const FaceInfo &face)
 Integration of local linear form.
 
virtual void face_residual (FEVector &face_vector1, FEVector &face_vector2, const FaceInfo &face1, const FaceInfo &face2)
 Integration of local linear form.
 
template<class BOX , class WORKER >
void integrate_1form (BOX &box, WORKER &local_forms, FEVector &dst, const FEVectors &src)
 The function integrating 1-forms as for instance used by residual().
 
void integrate_functional (LocalEstimate< dim > &local_forms, FEVectors &dst, const FEVectors &src)
 The function integrating functionals using local integrators.
 
void store_triangulation (Triangulation< dim > &new_tr)
 Function to copy a triangulation object for use after refinement.
 
void transfer_solution_to_coarse_mesh (Triangulation< dim > &tr_coarse, AppFrame::FEVector &coarse_solution, AppFrame::FEVector &refined_solution)
 Function to perform the transfer of a solution on a refined grid to the initial coarse grid.
 
void read_init_solution (AppFrame::FEVector &dst, bool &good_solution) const
 Function to transfer a solution on a refined grid to the initial coarse grid.
 
template<class INFO , typename TYPE >
void fill_local_data (const INFO &info, const std::vector< VectorSelector > &data_vector, std::vector< std::vector< std::vector< TYPE > > > &data) const
 Fill local data vector with values of the stored finite element function vectors.
 
virtual void grid_out (const std::string &basename)
 Output the grid used to solve the problem.
 
virtual void data_out (const std::string &basename, const FEVectors &src)
 
virtual void data_out (const std::string &basename, const FEVectors &src, const std::vector< std::string > &solution_names)
 
virtual void data_out (const std::string &basename, const FEVector &solution, const std::vector< std::string > &solution_names, const FEVector &postprocessing=FEVector(), const std::vector< std::string > &postprocessing_names=std::vector< std::string >())
 This function outputs the results of a computation.
 
- Public Member Functions inherited from AppFrame::ApplicationBase
 ApplicationBase (boost::shared_ptr< ApplicationData > ex_data=boost::shared_ptr< ApplicationData >())
 Constructor for an application.
 
 ApplicationBase (const ApplicationBase &other)
 Copy constructor.
 
virtual ~ApplicationBase ()
 Virtual destructor.
 
void print_parameters_to_file (ParameterHandler &param, const std::string &file_name, const ParameterHandler::OutputStyle &style)
 Print default parameters for the application to a file.
 
virtual void start_vector (FEVector &dst, std::string caller) const
 Initialize vector to problem size.
 
virtual void solve (FEVector &start, const FEVectors &rhs)
 Solve the system assembled with right hand side rhs and return the result in start.
 
virtual void Tsolve (FEVector &start, const FEVectors &rhs)
 Solve the dual system assembled with right hand side rhs and return the result in start.
 
virtual void grid_out (const std::string &filename) const
 Write the mesh in the format specified by the ParameterHandler.
 
virtual void data_out (const std::string &filename, const FEVectors &src, const std::vector< std::string >)
 
boost::shared_ptr
< ApplicationData
get_data ()
 Get access to the protected variable data.
 
const boost::shared_ptr
< ApplicationData
get_data () const
 Get read-only access to the protected variable data.
 
virtual std::string id () const
 Return a unique identification string for this application.
 
virtual void notify (const Event &reason)
 Add a reason for assembling.
 
virtual void clear ()
 Reset the application class such that a call to initialize() is possible again and will produce the same result as if the object was fresh.
 
void clear_events ()
 Clear all notifications.
 

Protected Attributes

Table< 2, DoFTools::Coupling > cell_couplings
 Couplings through cell bilinear forms.
 
Table< 2, DoFTools::Coupling > flux_couplings
 Couplings through flux bilinear forms.
 
Table< 2, DoFTools::Coupling > mg_cell_couplings
 Same as cell_couplings, but for multilevel matrices.
 
Table< 2, DoFTools::Coupling > mg_flux_couplings
 Same as flux_couplings, but for multilevel matrices.
 
boost::shared_ptr< Quadrature
< dim > > 
quadrature_assemble_cell
 Quadrature rule for matrix assembling on cells.
 
boost::shared_ptr< Quadrature
< dim-1 > > 
quadrature_assemble_face
 Quadrature rule for matrix assembling on faces.
 
std::vector< boost::shared_ptr
< MatrixType > > 
matrices
 Storage for the actual matrices.
 
std::vector< boost::shared_ptr
< MGMatrixType > > 
mg_matrices
 Storage for multilevel matrices.
 
std::vector< boost::shared_ptr
< MGMatrixType > > 
mg_interface_matrices_up
 The matrices at the interface between two refinement levels, coupling coarse to fine.
 
std::vector< boost::shared_ptr
< MGMatrixType > > 
mg_interface_matrices_down
 The matrices at the interface between two refinement levels, coupling fine to coarse.
 
- Protected Attributes inherited from AppFrame::DoFApplication< dim >
std::map< unsigned int, double > boundary_constraints
 List of all nodes constrained by a strong boundary condition, together with a value to be assigned.
 
ConstraintMatrix hanging_node_constraints
 Constraint Matrix object.
 
boost::shared_ptr< Mapping< dim > > mapping
 The mapping used for the transformation of mesh cells.
 
unsigned int mapping_degree
 Degree used for polynomial mapping of arbitrary order.
 
boost::shared_ptr
< Triangulation< dim > > 
tr
 Pointer to the Triangulation object.
 
boost::shared_ptr
< FiniteElement< dim > > 
element
 The finite element used in dof.
 
boost::shared_ptr< DoFHandler
< dim > > 
dof
 Pointer to the DoFHandler object.
 
boost::shared_ptr
< MGDoFHandler< dim > > 
mg_dof
 Pointer to the MGDoFHandler object.
 
MeshWorker::BlockInfo block_info
 
Vector< float > cell_errors
 The result of error estimation by cell.
 
Vector< float > face_errors
 The result of error estimation by face.
 
GridOut g_out
 The object for writing grids.
 
DataOut< dim, DoFHandler< dim > > d_out
 The object for writing data.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
data_interpretation
 Vector for adjusting output to vector data.
 
std::string refinement
 Refinement parameter from parameter file.
 
unsigned int initial_refinement
 Initial refinement from parameter file.
 
double refinement_threshold
 Refinement threshold for adaptive method from parameter file.
 
double coarsening_threshold
 Coarsening threshold for adaptive method from parameter file.
 
bool sort_cuthill
 Flag for sorting with Cuthill McKee algorithm.
 
Point< dimsort_direction
 Direction for downstream sorting.
 
std::vector< FEVector * > transfer_vectors
 List of vector names to be transfered from one grid to the next.
 
Quadrature< dimquadrature_residual_cell
 Quadrature rule for residual computation on cells.
 
Quadrature< dim-1 > quadrature_residual_bdry
 Quadrature rule for residual computation on boundary faces.
 
Quadrature< dim-1 > quadrature_residual_face
 Quadrature rule for residual computation on faces.
 
bool boundary_fluxes
 Extend the integration loops in assemble() and residual() also to boundary faces.
 
bool interior_fluxes
 Extend the integration loops in assemble() and residual() also to interior faces.
 
- Protected Attributes inherited from AppFrame::ApplicationBase
boost::shared_ptr
< ApplicationData
data
 Object for auxiliary data.
 
Event notifications
 Accumulate reasons for reassembling here.
 

Private Attributes

BlockSparsityPattern sparsity
 Sparsity patterns.
 
MGLevelObject
< BlockSparsityPattern > 
mg_sparsity
 Multi-grid sparsity patterns.
 
MGLevelObject
< BlockSparsityPattern > 
mg_sparsity_interface
 Sparsity patterns for #mg_interface_matrices in multigrid.
 
MeshWorker::VectorSelector cell_matrix_values
 The set of all vector data to be used in cell_matrix().
 
MeshWorker::VectorSelector flux_matrix_values
 The set of all vector data to be used in bdry_matrix() and face_matrix().
 

Additional Inherited Members

- Public Attributes inherited from AppFrame::DoFApplication< dim >
unsigned int verbosity
 Controls verbosity of certain functions.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
solution_interpretations
 solution_interpretations identifies whether some solution_names are scalars or parts of a vector.
 
std::vector
< DataComponentInterpretation::DataComponentInterpretation > 
postprocessing_interpretations
 postprocessing_interpretations identifies whether some postprocessing_names are scalars or parts of a vector.
 
bool output_materials_and_levels
 output_materials_and_levels if true then visualized, otherwise suppressed.
 
Vector< double > output_materials
 Vector that will be used by data_out to store the material ids.
 
Vector< double > output_levels
 Vector that will be used by data_out to store the refinement levels at each cell.
 
- Protected Member Functions inherited from AppFrame::DoFApplication< dim >
void constrain_boundary (FEVector &v, bool homogeneous) const
 Apply either homogeneous or inhomogeneous boundary_constraints.
 

Detailed Description

template<int dim>
class AppFrame::MatrixApplication< dim >

Application handling matrices and assembling linear systems of equations.

This class inherits the information on Triangulation and DoFHandler from its base class DoFApplication. It adds information on the linear system, in particular the sparsity.

Furthermore, it contains a vector of matrix blocks. This vector should be filled with void matrices and suitable block coordinates by the initialization routine of a derived class. After that, remesh_matrices() will reinitialize these matrices with appropriate sparsity patterns after each remesh().

Finally, generic loops over mesh cells and faces are implemented in order to assemble() the system matrices and compute the residual(). These functions use the virtual local integration functions declared here and implemented by the derived actual application class.

Note
You probably miss functions for building the right hand side of a linear system. They are missing deliberately, since the linear case is very special and can be handled in the framework of nonlinear problems. So, if you have a linear problem, code the local matrices for residual() so they ignore the BlockVector src.

Deriving your own application

Writing your own application derived from MatrixApplication requires setting up some data either in the constructor of your application or in its _initialize() function. Furthermore, some virtual functions should be implemented. Here's the list of things to do:

Providing information to MatrixApplication

MatrixApplication automizes the process of computing residuals and matrices. In order to do this efficiently, it requires some data structures to be set up in the beginning.

Note
Those of the following data not depending on the parameter file should be set up in the constructor, while your function _initialize() can be used to set up those depending on these parameters.
  1. All applications must fill the fields cell_couplings, flux_couplings if face integrals enter the formulation. These two objects are essentiel for saving memory and computation time.

  2. If multilevel support is required, then also mg_cell_couplings and mg_flux_couplings must be filled as well.

  3. If the formulation does not include integrals over interior or boundary faces, the variables boundary_fluxes and interior_fluxes should be set to false to save computation time.

  4. You must add matrices that will be filled by assemble() using add_matrix(). If using multigrid, do the same with add_mg_matrix().

  5. Finally, for problems depending on a previous finite element solution (e.g. in Newton's method or a timestepping scheme), you must tell MatrixApplication the name of the vector to be filled into the CellInfo and FaceInfo objects. This is done using the functions ???.

Author
Guido Kanschat, 2006, 2007

Member Typedef Documentation

The assembler object.

template<int dim>
typedef MatrixBlock<SparseMatrix<double> > AppFrame::MatrixApplication< dim >::MatrixType

The type of matrices stored.

The multigrid assembler object.

template<int dim>
typedef MatrixBlock<MGLevelObject<SparseMatrix<double> > > AppFrame::MatrixApplication< dim >::MGMatrixType

The type of multigrid matrices stored.

Constructor & Destructor Documentation

template<int dim>
AppFrame::MatrixApplication< dim >::MatrixApplication ( boost::shared_ptr< ApplicationData data = boost::shared_ptr< ApplicationData >(),
bool  multigrid = false 
)

Constructor for an object, owning its own mesh and dof handler.

If data is null, a new handler for application data is generated.

see constructor of DoFApplication.

template<int dim>
AppFrame::MatrixApplication< dim >::MatrixApplication ( bool  multigrid)

Constructor for an object, owning its own mesh and dof handler and having its own data.

template<int dim>
AppFrame::MatrixApplication< dim >::MatrixApplication ( DoFApplication< dim > &  ,
bool  triangulation_only,
bool  multigrid = false 
)

Constructor for an object, borrowing mesh and dof handler from another object.

If triangulation_only is true, only the triangulation is borrowed.

see constructor of DoFApplication.

Member Function Documentation

template<int dim>
void AppFrame::MatrixApplication< dim >::_initialize ( ParameterHandler &  param)

Initialize data of this class.

Note
remesh_dofs() has to be called before this member function since remesh_dofs() initializes data that is used in this routine. In particular, dof and block_indices.
template<int dim>
void AppFrame::MatrixApplication< dim >::add_matrix ( unsigned int  row,
unsigned int  col 
)

Add a new matrix to the list of matrices.

In order to initialize the matrix correctly, MatrixApplication needs to know the row and column this matrix is in.

Consecutive calls of add_matrix() generate a sequence of matrices in the same order/

template<int dim>
void AppFrame::MatrixApplication< dim >::add_mg_matrix ( unsigned int  row,
unsigned int  col 
)

Add a new multigrid matrix to mg_matrices.

See add_matrix() for details.

template<int dim>
void AppFrame::MatrixApplication< dim >::add_vector_for_cell_matrix ( std::string  name,
unsigned int  block,
unsigned int  components,
bool  values,
bool  derivatives 
)

Add a named vector from data to be handed to the cell_matrix() function via CellInfo::values and CellInfo::derivatives.

template<int dim>
void AppFrame::MatrixApplication< dim >::add_vector_for_flux_matrix ( std::string  name,
unsigned int  block,
unsigned int  components,
bool  values,
bool  derivatives 
)

Add a named vector from data to be handed to the bdry_matrix() and face_matrix() functions via FaceInfo::values and FaceInfo::derivatives.

template<int dim>
void AppFrame::MatrixApplication< dim >::assemble ( const FEVectors )

Loop over all cells and assemble the system matrices by using the local matrix integration functions cell_matrix(), bdry_matrix() and face_matrix() provided by the derived class.

Referenced by AppFrame::MatrixApplication< dim >::LocalMatrixIntegrator< ASSEMBLER >::bdry(), AppFrame::MatrixApplication< dim >::LocalMatrixIntegrator< ASSEMBLER >::cell(), and AppFrame::MatrixApplication< dim >::LocalMatrixIntegrator< ASSEMBLER >::face().

Here is the caller graph for this function:

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::bdry_matrix ( MatrixVector face_matrices,
const typename DoFApplication< dim >::FaceInfo face 
)
virtual

Integration of local bilinear form.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::cell_matrix ( MatrixVector cell_matrices,
const typename DoFApplication< dim >::CellInfo cell 
)
virtual

Integration of local bilinear form.

Here we loop over the quadrature points and over degrees of freedom in order to compute the matrix for the cell.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::declare_parameters ( ParameterHandler &  param)
virtual

Declare parameters related to matrices.

Note
The function declare_parameters of derived applications should always call the one of its base class.

Reimplemented from AppFrame::DoFApplication< dim >.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::face_matrix ( MatrixVector matrices11,
MatrixVector matrices12,
MatrixVector matrices21,
MatrixVector matrices22,
const typename DoFApplication< dim >::FaceInfo face1,
const typename DoFApplication< dim >::FaceInfo face2 
)
virtual

Integration of local bilinear form.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::initialize ( ParameterHandler &  param)
virtual

The virtual initialize function to be replaced by derived applications.

Reimplemented from AppFrame::DoFApplication< dim >.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::mg_bdry_matrix ( MatrixVector face_matrices,
const typename DoFApplication< dim >::FaceInfo face 
)
virtual

Integration of local bilinear form.

Note
Default implementation calls bdry_matrix.
template<int dim>
virtual void AppFrame::MatrixApplication< dim >::mg_cell_matrix ( MatrixVector cell_matrices,
const typename DoFApplication< dim >::CellInfo cell 
)
virtual

Integration of local bilinear form.

Here we loop over the quadrature points and over degrees of freedom in order to compute the matrix for the cell.

Note
Default implementation calls cell_matrix.
template<int dim>
virtual void AppFrame::MatrixApplication< dim >::mg_face_matrix ( MatrixVector matrices11,
MatrixVector matrices12,
MatrixVector matrices21,
MatrixVector matrices22,
const typename DoFApplication< dim >::FaceInfo face1,
const typename DoFApplication< dim >::FaceInfo face2 
)
virtual

Integration of local bilinear form.

Note
Default implementation calls face_matrix.
template<int dim>
virtual void AppFrame::MatrixApplication< dim >::post_cell_assemble ( )
virtual

A call back function for assemble(), called after all cell matrices have been entered, but before the face matrices are computed.

template<int dim>
virtual void AppFrame::MatrixApplication< dim >::remesh ( )
virtual

Refine grid accordingly to the Refinement options under "Grid Generation" in the parameter file.

Reimplemented from AppFrame::DoFApplication< dim >.

template<int dim>
void AppFrame::MatrixApplication< dim >::remesh_matrices ( )

Initialize sparsity patterns and matrices for the new mesh.

Member Data Documentation

template<int dim>
Table<2, DoFTools::Coupling> AppFrame::MatrixApplication< dim >::cell_couplings
protected

Couplings through cell bilinear forms.

cell_coupling allows to specify which variables couple with which equations. This is used by DoFTools to generate a sparsity pattern that does not contain elements unless specified in the cell_coupings, therefore, reducing the amount of memory needed.

Note
cell_coupings is a two dimensional table (i.e. a matrix) of DoFTools::Coupling objects. DoFTools::Coulings takes the values: always (components couple), zero (no coupling) and nonzero

Example: For the 2D Stokes equations

\begin{eqnarray*} \frac{d^2u}{dx^2} +\frac{d^2u}{dy^2} + \frac{dp}{dx} = 0 \\ \frac{d^2v}{dx^2} +\frac{d^2v}{dy^2} + \frac{dp}{dy} = 0 \\ \frac{du}{dx} +\frac{dv}{dy} = 0 \end{eqnarray*}

we need to initialize cell_coupling to:

* cell_coupling.reinit(3,3);
* cell_coupling(0,0) = DoFTools::always;
* cell_couplings(0,2) = DoFTools::always;
* cell_coupling(1,1) = DoFTools::always;
* ...
* 
template<int dim>
MeshWorker::VectorSelector AppFrame::MatrixApplication< dim >::cell_matrix_values
private

The set of all vector data to be used in cell_matrix().

template<int dim>
Table<2, DoFTools::Coupling> AppFrame::MatrixApplication< dim >::flux_couplings
protected

Couplings through flux bilinear forms.

template<int dim>
MeshWorker::VectorSelector AppFrame::MatrixApplication< dim >::flux_matrix_values
private

The set of all vector data to be used in bdry_matrix() and face_matrix().

template<int dim>
std::vector<boost::shared_ptr<MatrixType> > AppFrame::MatrixApplication< dim >::matrices
protected

Storage for the actual matrices.

This vector will be initialized by the concrete application, knowing how many matrices in each component are required.

It is suggested to use add_matrix() instead of initializing this array by hand.

template<int dim>
Table<2, DoFTools::Coupling> AppFrame::MatrixApplication< dim >::mg_cell_couplings
protected

Same as cell_couplings, but for multilevel matrices.

template<int dim>
Table<2, DoFTools::Coupling> AppFrame::MatrixApplication< dim >::mg_flux_couplings
protected

Same as flux_couplings, but for multilevel matrices.

template<int dim>
std::vector<boost::shared_ptr<MGMatrixType> > AppFrame::MatrixApplication< dim >::mg_interface_matrices_down
protected

The matrices at the interface between two refinement levels, coupling fine to coarse.

template<int dim>
std::vector<boost::shared_ptr<MGMatrixType> > AppFrame::MatrixApplication< dim >::mg_interface_matrices_up
protected

The matrices at the interface between two refinement levels, coupling coarse to fine.

template<int dim>
std::vector<boost::shared_ptr<MGMatrixType> > AppFrame::MatrixApplication< dim >::mg_matrices
protected

Storage for multilevel matrices.

template<int dim>
MGLevelObject<BlockSparsityPattern> AppFrame::MatrixApplication< dim >::mg_sparsity
private

Multi-grid sparsity patterns.

template<int dim>
MGLevelObject<BlockSparsityPattern> AppFrame::MatrixApplication< dim >::mg_sparsity_interface
private

Sparsity patterns for #mg_interface_matrices in multigrid.

template<int dim>
boost::shared_ptr<Quadrature<dim> > AppFrame::MatrixApplication< dim >::quadrature_assemble_cell
protected

Quadrature rule for matrix assembling on cells.

template<int dim>
boost::shared_ptr<Quadrature<dim-1> > AppFrame::MatrixApplication< dim >::quadrature_assemble_face
protected

Quadrature rule for matrix assembling on faces.

template<int dim>
BlockSparsityPattern AppFrame::MatrixApplication< dim >::sparsity
private

Sparsity patterns.


The documentation for this class was generated from the following file: