Galeri  Development
 All Classes Files Functions Variables Pages
Public Member Functions | List of all members
Galeri::FiniteElements::AbstractProblem Class Referenceabstract

Abstract interface to define linear problems. More...

#include <Galeri_AbstractProblem.h>

Inheritance diagram for Galeri::FiniteElements::AbstractProblem:
Inheritance graph
[legend]

Public Member Functions

virtual ~AbstractProblem ()
 Destructor.
 
virtual Epetra_RowMatrix & A ()=0
 Returns a reference to the linear system matrix.
 
virtual Epetra_MultiVector & RHS ()=0
 Returns a reference to the multi-vector of right-hand side.
 
virtual Epetra_MultiVector & LHS ()=0
 Returns a reference to the multi-vector of starting solution.
 
virtual const AbstractGridGrid () const =0
 Returns a reference to the grid object.
 
virtual const AbstractVariationalVariational () const =0
 Returns a reference to the variational object.
 
virtual void Compute ()=0
 Computes the linear system matrix, LHS and RHS.
 
virtual void ComputeNorms (Epetra_MultiVector &RowMatrixField, int(*ExactSolution)(double, double, double, double *), const bool verbose=true, double *SolutionNorm=0, double *ExactNorm=0, double *DiffNorm=0)=0
 Computes the norm of computed solution, exact solution, and error. More...
 

Detailed Description

Abstract interface to define linear problems.

AbstractProblem defines a set of abstract interfaces, used to construct the linear system corresponding to the finite element discretization of a scalar PDE problem. A concrete implementation will require an AbstractGrid and an AbstractVariational object; the former is used to query for the grid elements, the latter to integrate the variational form over such elements. The role of AbstractProblem is to take the elemental matrices, given by AbstractVariational, and insert them into the global, distributed matrix (whose RowMatrixRowMap() is given by Grid().RowMap()).

Author
Marzio Sala, SNL 9214.
Date
Last updated on Apr-05.

Member Function Documentation

virtual void Galeri::FiniteElements::AbstractProblem::ComputeNorms ( Epetra_MultiVector &  RowMatrixField,
int(*)(double, double, double, double *)  ExactSolution,
const bool  verbose = true,
double *  SolutionNorm = 0,
double *  ExactNorm = 0,
double *  DiffNorm = 0 
)
pure virtual

Computes the norm of computed solution, exact solution, and error.

Parameters
RowMatrixField- (In) Multi-vector defined on Grid().RowMap() which contains the numerical solution.
ExactSolution- (In) Function defined as in the following example: ExactSolution(double x, double y, double x, double* sol) will contain the value of the solution in sol[0], the x-derivative in sol[1], the y-derivative in sol[2], and the z-derivative in sol[2].
verbose- (In) If true, prints out the results.
SolutionNorm- (Out) a double array of size 3, which will contain the L2 norm, the semi-H1 norm and the H1-norm of the numerical solution.
ExactNorm- (Out) a double array of size 3, which will contain the L2 norm, the semi-H1 norm and the H1-norm of the exact solution.
DiffNorm- (Out) a double array of size 3, which will contain the L2 norm, the semi-H1 norm and the H1-norm of the error.

Implemented in Galeri::FiniteElements::LinearProblem.


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