BEM++
2.0
|
Function defined on a grid. More...
#include </home/wojtek/Projects/BEM/bempp-sven/bempp/lib/assembly/grid_function.hpp>
Public Types | |
enum | DataType { COEFFICIENTS, PROJECTIONS } |
enum | ConstructionMode { APPROXIMATE, INTERPOLATE } |
typedef Fiber::ScalarTraits < ResultType >::RealType | CoordinateType |
typedef Fiber::ScalarTraits < ResultType >::RealType | MagnitudeType |
Public Member Functions | |
GridFunction () | |
Constructor. More... | |
GridFunction (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const arma::Col< ResultType > &coefficients) | |
GridFunction (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const shared_ptr< const Space< BasisFunctionType > > &dualSpace, const arma::Col< ResultType > &projections) | |
GridFunction (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const shared_ptr< const Space< BasisFunctionType > > &dualSpace, const Function< ResultType > &function, ConstructionMode mode=APPROXIMATE) | |
Constructor. More... | |
BEMPP_DEPRECATED | GridFunction (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const shared_ptr< const Space< BasisFunctionType > > &dualSpace, const arma::Col< ResultType > &data, DataType dataType) |
Constructor. More... | |
BEMPP_DEPRECATED | GridFunction (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const shared_ptr< const Space< BasisFunctionType > > &dualSpace, const arma::Col< ResultType > &coefficients, const arma::Col< ResultType > &projections) |
Constructor. More... | |
bool | isInitialized () const |
Return whether this function has been properly initialized. | |
bool | wasInitializedFromCoefficients () const |
Returns true if the GridFunction was initialized with a list of coefficients and false if it was initialized with a list of projections. More... | |
shared_ptr< const Grid > | grid () const |
Grid on which this function is defined. More... | |
shared_ptr< const Space < BasisFunctionType > > | space () const |
Space in which this function is expanded. | |
shared_ptr< const Space < BasisFunctionType > > | dualSpace () const |
Space dual to the space in which this function is expanded. More... | |
shared_ptr< const Context < BasisFunctionType, ResultType > > | context () const |
Assembly context used to retrieve the strategy for evaluating any necessary integrals. | |
int | componentCount () const |
Number of components of this function. More... | |
const arma::Col< ResultType > & | coefficients () const |
Vector of expansion coefficients of this function in the basis of its expansion space. More... | |
arma::Col< ResultType > | projections (const shared_ptr< const Space< BasisFunctionType > > &dualSpace_) const |
Vector of scalar products of this function with the basis functions of dualSpace . More... | |
void | setCoefficients (const arma::Col< ResultType > &coeffs) |
Reset the expansion coefficients of this function in the basis of its primal space. More... | |
void | setProjections (const shared_ptr< const Space< BasisFunctionType > > &dualSpace_, const arma::Col< ResultType > &projects) |
Reinitialize the function by specifying the vector of its scalar products with the basis functions of dualSpace . More... | |
MagnitudeType | L2Norm () const |
Return the ![]() | |
const Fiber::Shapeset < BasisFunctionType > & | shapeset (const Entity< 0 > &element) const |
Return the shapeset associated with the given element. More... | |
const BEMPP_DEPRECATED Fiber::Basis < BasisFunctionType > & | basis (const Entity< 0 > &element) const |
Return the shapeset associated with the given element. More... | |
void | getLocalCoefficients (const Entity< 0 > &element, std::vector< ResultType > &coeffs) const |
Retrieve the expansion coefficients of this function on a single element. More... | |
void | evaluateAtSpecialPoints (VtkWriter::DataType dataType, arma::Mat< ResultType > &result_) const |
Evaluate function at either vertices or barycentres. More... | |
void | evaluateAtSpecialPoints (VtkWriter::DataType dataType, arma::Mat< CoordinateType > &points, arma::Mat< ResultType > &values) const |
void | evaluate (const Entity< 0 > &element, const arma::Mat< CoordinateType > &local, arma::Mat< ResultType > &values) const |
Evaluate function at specific points lying on a given element. More... | |
BEMPP_DEPRECATED arma::Col < ResultType > | projections (const Space< BasisFunctionType > &dualSpace_) const |
Vector of scalar products of this function with the basis functions of dualSpace . More... | |
BEMPP_DEPRECATED arma::Col < ResultType > | projections () const |
Vector of scalar products of this function with the basis functions of its dual space. More... | |
BEMPP_DEPRECATED void | setProjections (const Space< BasisFunctionType > &dualSpace_, const arma::Col< ResultType > &projects) |
Reinitialize the function by specifying the vector of its scalar products with the basis functions of dualSpace . More... | |
BEMPP_DEPRECATED void | setProjections (const arma::Col< ResultType > &projects) |
Reset the vector of scalar products of this function with the basis functions of its dual space. More... | |
BEMPP_DEPRECATED void | exportToVtk (VtkWriter::DataType dataType, const char *dataLabel, const char *fileNamesBase, const char *filesPath=0, VtkWriter::OutputType type=VtkWriter::ASCII) const |
Export this function to a VTK file. More... | |
BEMPP_DEPRECATED void | exportToGmsh (const char *dataLabel, const char *fileName) const |
Export this function to a Gmsh (.msh) file. More... | |
Private Member Functions | |
void | initializeFromCoefficients (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const arma::Col< ResultType > &coefficients) |
void | initializeFromProjections (const shared_ptr< const Context< BasisFunctionType, ResultType > > &context, const shared_ptr< const Space< BasisFunctionType > > &space, const shared_ptr< const Space< BasisFunctionType > > &dualSpace, const arma::Col< ResultType > &projections) |
void | updateProjectionsFromCoefficients (const shared_ptr< const Space< BasisFunctionType > > &dualSpace_) const |
void | updateCoefficientsFromProjections () const |
Private Attributes | |
shared_ptr< const Context < BasisFunctionType, ResultType > > | m_context |
shared_ptr< const Space < BasisFunctionType > > | m_space |
shared_ptr< const Space < BasisFunctionType > > | m_dualSpace |
shared_ptr< const arma::Col < ResultType > > | m_coefficients |
shared_ptr< const arma::Col < ResultType > > | m_projections |
bool | m_wasInitializedFromCoefficients |
Related Functions | |
(Note that these are not member functions.) | |
template<typename BasisFunctionType , typename ResultType > | |
GridFunction < BasisFunctionType, ResultType > | operator+ (const GridFunction< BasisFunctionType, ResultType > &g) |
Return a copy of the passed function g . | |
template<typename BasisFunctionType , typename ResultType > | |
GridFunction < BasisFunctionType, ResultType > | operator- (const GridFunction< BasisFunctionType, ResultType > &g) |
Return the grid function representing the function g multipled by -1. More... | |
template<typename BasisFunctionType , typename ResultType > | |
GridFunction < BasisFunctionType, ResultType > | operator+ (const GridFunction< BasisFunctionType, ResultType > &g1, const GridFunction< BasisFunctionType, ResultType > &g2) |
Return a grid function representing the sum of the operands. More... | |
template<typename BasisFunctionType , typename ResultType > | |
GridFunction < BasisFunctionType, ResultType > | operator- (const GridFunction< BasisFunctionType, ResultType > &g1, const GridFunction< BasisFunctionType, ResultType > &g2) |
Return a grid function representing the difference of the operands. More... | |
template<typename BasisFunctionType , typename ResultType , typename ScalarType > | |
GridFunction < BasisFunctionType, ResultType > | operator* (const GridFunction< BasisFunctionType, ResultType > &g, const ScalarType &scalar) |
Return the grid function representing the function g multiplied by scalar . More... | |
template<typename BasisFunctionType , typename ResultType , typename ScalarType > | |
boost::enable_if< typename boost::mpl::has_key < boost::mpl::set< float, double, std::complex< float > , std::complex< double > >, ScalarType >, GridFunction < BasisFunctionType, ResultType > >::type | operator* (const ScalarType &scalar, const GridFunction< BasisFunctionType, ResultType > &g) |
Return the grid function representing the function g multiplied by scalar . More... | |
template<typename BasisFunctionType , typename ResultType , typename ScalarType > | |
GridFunction < BasisFunctionType, ResultType > | operator/ (const GridFunction< BasisFunctionType, ResultType > &g1, const ScalarType &scalar) |
Return the grid function representing the function g divided by scalar . More... | |
template<typename BasisFunctionType , typename ResultType > | |
void | exportToVtk (const GridFunction< BasisFunctionType, ResultType > &gridFunction, VtkWriter::DataType dataType, const char *dataLabel, const char *fileNamesBase, const char *filesPath=0, VtkWriter::OutputType type=VtkWriter::ASCII) |
Export this function to a VTK file. More... | |
template<typename BasisFunctionType , typename ResultType > | |
void | exportToGmsh (const GridFunction< BasisFunctionType, ResultType > &gridFunction, const char *dataLabel, const char *fileName) |
Export this function to a Gmsh (.msh) file. More... | |
template<typename BasisFunctionType , typename ResultType > | |
ScalarTraits < BasisFunctionType > ::RealType | L2NormOfDifference (const GridFunction< BasisFunctionType, ResultType > &gridFunction, const Fiber::Function< ResultType > &refFunction, const Fiber::QuadratureStrategy< BasisFunctionType, ResultType, GeometryFactory > &quadStrategy, const EvaluationOptions &options=EvaluationOptions()) |
Calculate the ![]() | |
template<typename BasisFunctionType , typename ResultType > | |
void | estimateL2Error (const GridFunction< BasisFunctionType, ResultType > &gridFunction, const Fiber::Function< ResultType > &refFunction, const Fiber::QuadratureStrategy< BasisFunctionType, ResultType, GeometryFactory > &quadStrategy, const EvaluationOptions &options, typename ScalarTraits< BasisFunctionType >::RealType &absError, typename ScalarTraits< BasisFunctionType >::RealType &relError) |
Calculate the absolute and relative ![]() | |
template<typename BasisFunctionType , typename ResultType > | |
void | estimateL2Error (const GridFunction< BasisFunctionType, ResultType > &gridFunction, const Fiber::Function< ResultType > &refFunction, const Fiber::QuadratureStrategy< BasisFunctionType, ResultType, GeometryFactory > &quadStrategy, typename ScalarTraits< BasisFunctionType >::RealType &absError, typename ScalarTraits< BasisFunctionType >::RealType &relError) |
Calculate the absolute and relative ![]() | |
Function defined on a grid.
This class represents a function defined on a grid and expanded in a particular function space.
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | ) |
Constructor.
Construct an uninitialized grid function. The only way to initialize it later is using the assignment operator.
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | const shared_ptr< const Context< BasisFunctionType, ResultType > > & | context, |
const shared_ptr< const Space< BasisFunctionType > > & | space, | ||
const arma::Col< ResultType > & | coefficients | ||
) |
Constructor.
[in] | context | Assembly context from which a quadrature strategy can be retrieved. |
[in] | space | Function space to expand the grid function in. |
[in] | coefficients | Vector of length space.globalDofCount() containing the expansion coefficients of the grid function in the space space . |
This constructor builds a grid function from its coefficients in a function space.
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | const shared_ptr< const Context< BasisFunctionType, ResultType > > & | context, |
const shared_ptr< const Space< BasisFunctionType > > & | space, | ||
const shared_ptr< const Space< BasisFunctionType > > & | dualSpace, | ||
const arma::Col< ResultType > & | projections | ||
) |
Constructor.
[in] | context | Assembly context from which a quadrature strategy can be retrieved. |
[in] | space | Function space to expand the grid function in. |
[in] | dualSpace | Function space dual to space . |
[in] | projections | Vector of length dualSpace.globalDofCount() containing the scalar products of the grid function and the basis functions of the space dualSpace . |
This constructor builds a grid function expanded in the basis of the space space
from the projections of this function on the basis functions of another space dualSpace
. Both spaces must be defined on the same grid.
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | const shared_ptr< const Context< BasisFunctionType, ResultType > > & | context, |
const shared_ptr< const Space< BasisFunctionType > > & | space, | ||
const shared_ptr< const Space< BasisFunctionType > > & | dualSpace, | ||
const Function< ResultType > & | function, | ||
ConstructionMode | mode = APPROXIMATE |
||
) |
Constructor.
[in] | context | Assembly context from which a quadrature strategy can be retrieved. |
[in] | space | Function space to expand the grid function in. |
[in] | dualSpace | Function space dual to space . |
[in] | function | Function object whose values on space.grid() will be used to construct the new grid function. |
[in] | mode | Determines the way the expansion coefficients of the function are determined. Can be set to #APPROXIMATE or #INTERPOLATE. These options are discussed below. |
This constructor builds a grid function belonging to the space space
and approximating the function defined by the object
function
.
If mode
is set to #APPROXIMATE, the approximate coefficients of
function
in the basis of
space
are determined by solving the equation
in the least-squares sense, with denoting the set of basis functions of
dualSpace
. The two spaces must be defined on the same grid.
Otherwise (if mode
is set to #INTERPOLATE), the approximate coefficients are determined by evaluating
function
at the interpolation points of the global degrees of freedom of space
(see Space::getGlobalDofInterpolationPoints() for more details).
References Bempp::GridFunction< BasisFunctionType, ResultType >::setCoefficients(), and Bempp::GridFunction< BasisFunctionType, ResultType >::setProjections().
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | const shared_ptr< const Context< BasisFunctionType, ResultType > > & | context, |
const shared_ptr< const Space< BasisFunctionType > > & | space, | ||
const shared_ptr< const Space< BasisFunctionType > > & | dualSpace, | ||
const arma::Col< ResultType > & | data, | ||
DataType | dataType | ||
) |
Constructor.
[in] | context | Assembly context from which a quadrature strategy can be retrieved. |
[in] | space | Function space to expand the grid function in. |
[in] | dualSpace | Function space dual to space . |
[in] | data | If dataType == COEFFICIENTS , the vector data should have length space.globalDofCount() and contain the expansion coefficients of the grid function in the space space . Otherwise, if dataType == PROJECTIONS , data should have length dualSpace.globalDofCount() contain the scalar products of the grid function and the basis functions of the space dualSpace . |
[in] | dataType | Interpretation of the vector passed via the argument data . |
This constructor builds a grid function from either its coefficients in a function space or from its projections on the basis functions of another (dual) function space. If the other type of data turns out to be needed later (for example, the projections vector if coefficients were supplied in the constructor), it is calculated automatically. The Context object given in the constructor is then used to determine the strategy for evaluating any necessary integrals.
space
and dualSpace
must be defined on the same grid.
Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction | ( | const shared_ptr< const Context< BasisFunctionType, ResultType > > & | context, |
const shared_ptr< const Space< BasisFunctionType > > & | space, | ||
const shared_ptr< const Space< BasisFunctionType > > & | dualSpace, | ||
const arma::Col< ResultType > & | coefficients, | ||
const arma::Col< ResultType > & | projections | ||
) |
Constructor.
[in] | context | Assembly context from which a quadrature strategy can be retrieved. |
[in] | space | Function space to expand the grid function in. |
[in] | dualSpace | Function space dual to space . |
[in] | coefficients | Vector of length space.globalDofCount() containing the expansion coefficients of the grid function in the space space . |
[in] | projections | Vector of length dualSpace.globalDofCount() containing the scalar products of the grid function and the basis functions of the space dualSpace . |
space
and dualSpace
must be defined on the same grid.
const Fiber::Basis< BasisFunctionType > & Bempp::GridFunction< BasisFunctionType, ResultType >::basis | ( | const Entity< 0 > & | element | ) | const |
Return the shapeset associated with the given element.
const arma::Col< ResultType > & Bempp::GridFunction< BasisFunctionType, ResultType >::coefficients | ( | ) | const |
Vector of expansion coefficients of this function in the basis of its expansion space.
An exception is thrown if this function is called on an uninitialized GridFunction object (one constructed with the default constructor).
Referenced by Bempp::AssembledPotentialOperator< BasisFunctionType, ResultType >::apply(), Bempp::BoundaryOperator< BasisFunctionType, ResultType >::apply(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator*(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator+(), and Bempp::GridFunction< BasisFunctionType, ResultType >::operator-().
int Bempp::GridFunction< BasisFunctionType, ResultType >::componentCount | ( | ) | const |
Number of components of this function.
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::exportToGmsh().
shared_ptr< const Space< BasisFunctionType > > Bempp::GridFunction< BasisFunctionType, ResultType >::dualSpace | ( | ) | const |
Space dual to the space in which this function is expanded.
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::operator*(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator+(), and Bempp::GridFunction< BasisFunctionType, ResultType >::operator-().
void Bempp::GridFunction< BasisFunctionType, ResultType >::evaluate | ( | const Entity< 0 > & | element, |
const arma::Mat< CoordinateType > & | local, | ||
arma::Mat< ResultType > & | values | ||
) | const |
Evaluate function at specific points lying on a given element.
[in] | element | An element belonging to the grid on which function is defined. |
[in] | local | A 2D array whose (i,j)th element is the ith coordinate of the jth point at which the grid function should be evaluated, in the local coordinate system of element . |
[in] | values | A 2D array whose (i,j)th element is the ith component of function evaluated at the jth point. |
References Fiber::CollectionOfShapesetTransformations< CoordinateType_ >::addDependencies(), Fiber::Shapeset< BasisFunctionType >::evaluate(), Fiber::CollectionOfShapesetTransformations< CoordinateType_ >::evaluate(), Bempp::Entity< 0 >::geometry(), Fiber::CollectionOfShapesetTransformations< CoordinateType_ >::resultDimension(), Fiber::Shapeset< BasisFunctionType >::size(), and Fiber::CollectionOfShapesetTransformations< CoordinateType_ >::transformationCount().
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::exportToGmsh().
void Bempp::GridFunction< BasisFunctionType, ResultType >::evaluateAtSpecialPoints | ( | VtkWriter::DataType | dataType, |
arma::Mat< ResultType > & | result_ | ||
) | const |
Evaluate function at either vertices or barycentres.
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::exportToVtk().
void Bempp::GridFunction< BasisFunctionType, ResultType >::exportToGmsh | ( | const char * | dataLabel, |
const char * | fileName | ||
) | const |
Export this function to a Gmsh (.msh) file.
[in] | dataLabel | Label used to identify the function in the Gmsh file. |
[in] | fileName | Name of the file to be written. |
.msh file generated by exportToGmsh cannot be imported back by BEM++ and used to create a Grid object.void Bempp::GridFunction< BasisFunctionType, ResultType >::exportToVtk | ( | VtkWriter::DataType | dataType, |
const char * | dataLabel, | ||
const char * | fileNamesBase, | ||
const char * | filesPath = 0 , |
||
VtkWriter::OutputType | type = VtkWriter::ASCII |
||
) | const |
Export this function to a VTK file.
[in] | dataType | Determines whether data are attaches to vertices or cells. |
[in] | dataLabel | Label used to identify the function in the VTK file. |
[in] | fileNamesBase | Base name of the output files. It should not contain any directory part or filename extensions. |
[in] | filesPath | Output directory. Can be set to NULL, in which case the files are output in the current directory. |
[in] | type | Output type (default: ASCII). See Dune reference manual for more details. |
void Bempp::GridFunction< BasisFunctionType, ResultType >::getLocalCoefficients | ( | const Entity< 0 > & | element, |
std::vector< ResultType > & | coeffs | ||
) | const |
Retrieve the expansion coefficients of this function on a single element.
[in] | element | An element belonging to the grid space.grid() . |
[out] | coeffs | Vector of the expansion coefficients of this function corresponding to the basis functions of the primal space living on element element . |
shared_ptr< const Grid > Bempp::GridFunction< BasisFunctionType, ResultType >::grid | ( | ) | const |
Grid on which this function is defined.
Referenced by Bempp::ElementaryPotentialOperator< BasisFunctionType_, KernelType_, ResultType_ >::evaluateAtPoints(), and Bempp::ElementaryPotentialOperator< BasisFunctionType_, KernelType_, ResultType_ >::evaluateOnGrid().
GridFunction< BasisFunctionType, ResultType >::MagnitudeType Bempp::GridFunction< BasisFunctionType, ResultType >::L2Norm | ( | ) | const |
Return the -norm of the grid function.
References Bempp::NO_TRANSPOSE.
arma::Col< ResultType > Bempp::GridFunction< BasisFunctionType, ResultType >::projections | ( | const shared_ptr< const Space< BasisFunctionType > > & | dualSpace_ | ) | const |
Vector of scalar products of this function with the basis functions of dualSpace
.
dualSpace
must be defined on the same grid as the space in which the GridFunction is expanded.
An exception is thrown if this function is called on an uninitialized GridFunction object (one constructed with the default constructor).
Referenced by Bempp::BoundaryOperator< BasisFunctionType, ResultType >::apply(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator*(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator+(), and Bempp::GridFunction< BasisFunctionType, ResultType >::operator-().
arma::Col< ResultType > Bempp::GridFunction< BasisFunctionType, ResultType >::projections | ( | const Space< BasisFunctionType > & | dualSpace_ | ) | const |
Vector of scalar products of this function with the basis functions of dualSpace
.
dualSpace
must be defined on the same grid as the space in which the GridFunction is expanded.
An exception is thrown if this function is called on an uninitialized GridFunction object (one constructed with the default constructor).
arma::Col< ResultType > Bempp::GridFunction< BasisFunctionType, ResultType >::projections | ( | ) | const |
Vector of scalar products of this function with the basis functions of its dual space.
An exception is thrown if this function is called on an uninitialized GridFunction object (one constructed with the default constructor).
void Bempp::GridFunction< BasisFunctionType, ResultType >::setCoefficients | ( | const arma::Col< ResultType > & | coeffs | ) |
Reset the expansion coefficients of this function in the basis of its primal space.
As a side effect, any internally stored vector of the projections of this grid function on the basis functions of its dual space is marked as invalid and recalculated on the next call to projections().
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction().
void Bempp::GridFunction< BasisFunctionType, ResultType >::setProjections | ( | const shared_ptr< const Space< BasisFunctionType > > & | dualSpace_, |
const arma::Col< ResultType > & | projects | ||
) |
Reinitialize the function by specifying the vector of its scalar products with the basis functions of dualSpace
.
dualSpace
must be defined on the same grid as the space in which the GridFunction is expanded.
Referenced by Bempp::BoundaryOperator< BasisFunctionType, ResultType >::apply(), and Bempp::GridFunction< BasisFunctionType, ResultType >::GridFunction().
void Bempp::GridFunction< BasisFunctionType, ResultType >::setProjections | ( | const Space< BasisFunctionType > & | dualSpace_, |
const arma::Col< ResultType > & | projects | ||
) |
Reinitialize the function by specifying the vector of its scalar products with the basis functions of dualSpace
.
dualSpace
must be defined on the same grid as the space in which the GridFunction is expanded.
dualSpace_
is not destroyed as long as this GridFunction object is alive. void Bempp::GridFunction< BasisFunctionType, ResultType >::setProjections | ( | const arma::Col< ResultType > & | projects | ) |
Reset the vector of scalar products of this function with the basis functions of its dual space.
const Fiber::Shapeset< BasisFunctionType > & Bempp::GridFunction< BasisFunctionType, ResultType >::shapeset | ( | const Entity< 0 > & | element | ) | const |
Return the shapeset associated with the given element.
bool Bempp::GridFunction< BasisFunctionType, ResultType >::wasInitializedFromCoefficients | ( | ) | const |
Returns true
if the GridFunction was initialized with a list of coefficients and false
if it was initialized with a list of projections.
This means that the other quantity (projections or coefficients, respectively) is a derived quantity and may therefore be less accurate.
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::operator*(), Bempp::GridFunction< BasisFunctionType, ResultType >::operator+(), and Bempp::GridFunction< BasisFunctionType, ResultType >::operator-().
|
related |
Calculate the absolute and relative errors of a solution.
This function calculates the absolute and relative norms of the difference between a GridFunction (typically a numerical solution of an integral equation) and a Function (typically a known analytical solution.
The quadrature strategy quadStrategy
is used to evaluate any necessary integrals; the options
object controls the level of parallelism. The results are returned in the output arguments absError
and relError
.
|
related |
Calculate the absolute and relative errors of a solution.
This is an overloaded version, provided for convenience. It is equivalent to the six-parameter version with options
set to EvaluationOptions()
.
|
related |
Export this function to a Gmsh (.msh) file.
[in] | dataLabel | Label used to identify the function in the Gmsh file. |
[in] | fileName | Name of the file to be written. |
.msh file generated by exportToGmsh cannot be imported back by BEM++ and used to create a Grid object. References Bempp::GridFunction< BasisFunctionType, ResultType >::componentCount(), Bempp::GridView::entityIterator(), Bempp::GridFunction< BasisFunctionType, ResultType >::evaluate(), Bempp::Entity< 0 >::geometry(), Bempp::Geometry::local2global(), Bempp::GridFunction< BasisFunctionType, ResultType >::space(), and Bempp::Entity< 0 >::type().
|
related |
Export this function to a VTK file.
[in] | dataType | Determines whether data are attaches to vertices or cells. |
[in] | dataLabel | Label used to identify the function in the VTK file. |
[in] | fileNamesBase | Base name of the output files. It should not contain any directory part or filename extensions. |
[in] | filesPath | Output directory. Can be set to NULL, in which case the files are output in the current directory. |
[in] | type | Output type (default: ASCII). See Dune reference manual for more details. |
References Bempp::GridFunction< BasisFunctionType, ResultType >::evaluateAtSpecialPoints(), Bempp::GridFunction< BasisFunctionType, ResultType >::space(), and Bempp::GridView::vtkWriter().
|
related |
Calculate the -norm of the difference between a GridFunction and a Function.
This function can be used to estimate the absolute error of a numerical solution gridFunction
against a known analytical solution refFunction
.
The quadrature strategy quadStrategy
is used to evaluate any necessary integrals; the options
object controls the level of parallelism.
|
related |
Return the grid function representing the function g
multiplied by scalar
.
g
is uninitialized. References Bempp::GridFunction< BasisFunctionType, ResultType >::coefficients(), Bempp::GridFunction< BasisFunctionType, ResultType >::context(), Bempp::GridFunction< BasisFunctionType, ResultType >::dualSpace(), Bempp::GridFunction< BasisFunctionType, ResultType >::projections(), Bempp::GridFunction< BasisFunctionType, ResultType >::space(), and Bempp::GridFunction< BasisFunctionType, ResultType >::wasInitializedFromCoefficients().
|
related |
Return the grid function representing the function g
multiplied by scalar
.
g
is uninitialized.
|
related |
Return a grid function representing the sum of the operands.
References Bempp::GridFunction< BasisFunctionType, ResultType >::coefficients(), Bempp::GridFunction< BasisFunctionType, ResultType >::context(), Bempp::GridFunction< BasisFunctionType, ResultType >::dualSpace(), Bempp::GridFunction< BasisFunctionType, ResultType >::projections(), Bempp::GridFunction< BasisFunctionType, ResultType >::space(), and Bempp::GridFunction< BasisFunctionType, ResultType >::wasInitializedFromCoefficients().
|
related |
Return the grid function representing the function g
multipled by -1.
g
is uninitialized.
|
related |
Return a grid function representing the difference of the operands.
References Bempp::GridFunction< BasisFunctionType, ResultType >::coefficients(), Bempp::GridFunction< BasisFunctionType, ResultType >::context(), Bempp::GridFunction< BasisFunctionType, ResultType >::dualSpace(), Bempp::GridFunction< BasisFunctionType, ResultType >::projections(), Bempp::GridFunction< BasisFunctionType, ResultType >::space(), and Bempp::GridFunction< BasisFunctionType, ResultType >::wasInitializedFromCoefficients().
|
related |
Return the grid function representing the function g
divided by scalar
.
g
is uninitialized.