BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Public Member Functions | Protected Member Functions | Private Member Functions | Related Functions | List of all members
Bempp::DiscreteBoundaryOperator< ValueType > Class Template Referenceabstract

Discrete boundary operator. More...

#include </home/wojtek/Projects/BEM/bempp-sven/bempp/lib/assembly/discrete_boundary_operator.hpp>

Inheritance diagram for Bempp::DiscreteBoundaryOperator< ValueType >:
Bempp::AcaApproximateLuInverse< ValueType > Bempp::DiscreteAcaBoundaryOperator< ValueType > Bempp::DiscreteBlockedBoundaryOperator< ValueType > Bempp::DiscreteBoundaryOperatorComposition< ValueType > Bempp::DiscreteBoundaryOperatorSum< ValueType > Bempp::DiscreteDenseBoundaryOperator< ValueType > Bempp::DiscreteInverseSparseBoundaryOperator< ValueType > Bempp::DiscreteNullBoundaryOperator< ValueType > Bempp::DiscreteSparseBoundaryOperator< ValueType > Bempp::ScaledDiscreteBoundaryOperator< ValueType > Bempp::TransposedDiscreteBoundaryOperator< ValueType >

Public Member Functions

virtual ~DiscreteBoundaryOperator ()
 Destructor.
 
void apply (const Thyra::EOpTransp M_trans, const Thyra::MultiVectorBase< ValueType > &x_in, const Thyra::Ptr< Thyra::MultiVectorBase< ValueType > > &y_inout, const ValueType alpha, const ValueType beta) const
 Apply the linear operator to a multivector. More...
 
void apply (const TranspositionMode trans, const arma::Mat< ValueType > &x_in, arma::Mat< ValueType > &y_inout, const ValueType alpha, const ValueType beta) const
 Apply the operator to a vector. More...
 
virtual shared_ptr< const
DiscreteBoundaryOperator
asDiscreteAcaBoundaryOperator (double eps=-1, int maximumRank=-1, bool interleave=false) const
 Return a representation that can be cast to a DiscreteAcaBoundaryOperator. More...
 
virtual void dump () const
 Write a textual representation of the operator to standard output. More...
 
virtual arma::Mat< ValueType > asMatrix () const
 Matrix representation of the operator. More...
 
virtual unsigned int rowCount () const =0
 Number of rows of the operator.
 
virtual unsigned int columnCount () const =0
 Number of columns of the operator.
 
virtual void addBlock (const std::vector< int > &rows, const std::vector< int > &cols, const ValueType alpha, arma::Mat< ValueType > &block) const =0
 Add a subblock of this operator to a matrix. More...
 

Protected Member Functions

virtual void applyImpl (const Thyra::EOpTransp M_trans, const Thyra::MultiVectorBase< ValueType > &X_in, const Teuchos::Ptr< Thyra::MultiVectorBase< ValueType > > &Y_inout, const ValueType alpha, const ValueType beta) const
 

Private Member Functions

virtual void applyBuiltInImpl (const TranspositionMode trans, const arma::Col< ValueType > &x_in, arma::Col< ValueType > &y_inout, const ValueType alpha, const ValueType beta) const =0
 

Related Functions

(Note that these are not member functions.)

template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
operator+ (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Unary plus: return a copy of the argument.
 
template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
operator- (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the operand multiplied by -1. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
operator+ (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2)
 Return a shared pointer to a DiscreteBoundaryOperator representing the sum of the operands. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
operator- (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2)
 Return a shared pointer to a DiscreteBoundaryOperator representing the difference of the operands. More...
 
template<typename ValueType , typename ScalarType >
boost::enable_if< typename
boost::mpl::has_key
< boost::mpl::set< float,
double, std::complex< float >
, std::complex< double >
>, ScalarType >, shared_ptr
< DiscreteBoundaryOperator
< ValueType > > >::type 
operator* (ScalarType scalar, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op multiplied by scalar. More...
 
template<typename ValueType , typename ScalarType >
boost::enable_if< typename
boost::mpl::has_key
< boost::mpl::set< float,
double, std::complex< float >
, std::complex< double >
>, ScalarType >, shared_ptr
< DiscreteBoundaryOperator
< ValueType > > >::type 
operator* (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op, ScalarType scalar)
 Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op multiplied by scalar. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
operator* (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2)
 Return a shared pointer to a DiscreteBoundaryOperator representing a composition of operators, (*op1) * (*op2). More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
mul (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2)
 Return a shared pointer to a DiscreteBoundaryOperator representing a composition of operators, (*op1) * (*op2). More...
 
template<typename ValueType , typename ScalarType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
operator/ (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op, ScalarType scalar)
 Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op divided by scalar. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
transpose (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the transposed operator *op. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
conjugate (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the conjugated operator *op. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
conjugateTranspose (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the conjugate-transposed operator *op. More...
 
template<typename ValueType >
shared_ptr
< DiscreteBoundaryOperator
< ValueType > > 
transpose (TranspositionMode trans, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the transposed and/or conjugated (depending on the value of the argument trans) operator *op. More...
 
template<typename RealType >
shared_ptr
< DiscreteBoundaryOperator
< std::complex< RealType > > > 
complexify (const shared_ptr< const DiscreteBoundaryOperator< RealType > > &op)
 Return a shared pointer to a DiscreteBoundaryOperator representing the "complexified" operator *op. More...
 

Detailed Description

template<typename ValueType>
class Bempp::DiscreteBoundaryOperator< ValueType >

Discrete boundary operator.

This class represents a discretised boundary operator, i.e. a matrix $\mathsf{L}$ with entries

\[ \mathsf{L}_{mn} \equiv \int_S \phi^*_m(x) \,[L\, \psi_n](x)\, \mathrm{d}S(x), \]

where $L$ is a boundary operator, $\{\phi_m\}_{m=1}^{M}$ are the basis functions spanning a test space and defined on a surface $S$, whereas $\{\psi_n\}_{n=1}^{N}$ are the basis functions of a trial space. The way the matrix is stored can be arbitrary. Concrete subclasses of this class implement specific storage methods.

If BEM++ is compiled with Trilinos, this class is derived from Thyra::LinearOpDefaultBase<ValueType> and hence inherits all non-private member function of this class; see Trilinos documentation for full documentation of these functions.

If Trilinos is not available during compilation, a simple fallback interface is provided.

Note
The lifetime of DiscreteBoundaryOperator objects must managed by shared pointers.
Template Parameters
ValueTypeType used to represent entries of the operator matrix. It can take the following values: float, double, std::complex<float> and std::complex<double>.

Member Function Documentation

template<typename ValueType>
virtual void Bempp::DiscreteBoundaryOperator< ValueType >::addBlock ( const std::vector< int > &  rows,
const std::vector< int > &  cols,
const ValueType  alpha,
arma::Mat< ValueType > &  block 
) const
pure virtual

Add a subblock of this operator to a matrix.

Perform the operation block += alpha * L[rows, cols], where block is a matrix, alpha a scalar and L[rows, cols] a subblock of this discrete linear operator.

Parameters
[in]rowsVector of row indices.
[in]colsVector of column indices.
[in]alphaMultiplier.
[in,out]blockOn entry, matrix of size (rows.size(), cols.size()). On exit, each element (i, j) of this matrix will be augmented by the element (rows[i], cols[j]) of this operator, multiplied by alpha.

Row and column indices may be unsorted.

This method need not be supported by all subclasses. It is mainly intended for internal use during weak-form assembly.

Implemented in Bempp::DiscreteAcaBoundaryOperator< ValueType >, Bempp::DiscreteSparseBoundaryOperator< ValueType >, Bempp::DiscreteBlockedBoundaryOperator< ValueType >, Bempp::DiscreteInverseSparseBoundaryOperator< ValueType >, Bempp::AcaApproximateLuInverse< ValueType >, Bempp::DiscreteBoundaryOperatorSum< ValueType >, Bempp::ScaledDiscreteBoundaryOperator< ValueType >, Bempp::ComplexifiedDiscreteBoundaryOperator< RealType >, Bempp::TransposedDiscreteBoundaryOperator< ValueType >, Bempp::DiscreteBoundaryOperatorComposition< ValueType >, Bempp::DiscreteDenseBoundaryOperator< ValueType >, and Bempp::DiscreteNullBoundaryOperator< ValueType >.

template<typename ValueType>
void Bempp::DiscreteBoundaryOperator< ValueType >::apply ( const Thyra::EOpTransp  M_trans,
const Thyra::MultiVectorBase< ValueType > &  x_in,
const Thyra::Ptr< Thyra::MultiVectorBase< ValueType > > &  y_inout,
const ValueType  alpha,
const ValueType  beta 
) const

Apply the linear operator to a multivector.

Set the elements of the multivector y_inout to

y_inout := alpha * trans(L) * x_in + beta * y_inout,

where L is the linear operator represented by this object.

Parameters
[in]M_transDetermines whether what is applied is the "bare" operator, its transpose, conjugate or conjugate transpose.
[in]x_inThe right-hand-side multivector.
[in,out]y_inoutThe target multivector being transformed. When beta == 0.0, this multivector can have uninitialized elements.
[in]alphaScalar multiplying this operator.
[in]betaThe multiplier for the target multivector y_inout.

This overload of the apply() member function is only available if the library was compiled with Trilinos.

template<typename ValueType>
void Bempp::DiscreteBoundaryOperator< ValueType >::apply ( const TranspositionMode  trans,
const arma::Mat< ValueType > &  x_in,
arma::Mat< ValueType > &  y_inout,
const ValueType  alpha,
const ValueType  beta 
) const

Apply the operator to a vector.

Set the elements of the column vector y_inout to

y_inout := alpha * trans(L) * x_in + beta * y_inout,

where L is the linear operator represented by this object.

Parameters
[in]transDetermines whether what is applied is the "bare" operator, its transpose, conjugate or conjugate transpose.
[in]x_inThe right-hand-side multivector.
[in,out]y_inoutThe target multivector being transformed. When beta == 0.0, this multivector can have uninitialized elements.
[in]alphaScalar multiplying this operator.
[in]betaThe multiplier for the target multivector y_inout.

This overload is always available, even if the library was compiled without Trilinos.

References Bempp::CONJUGATE_TRANSPOSE, and Bempp::TRANSPOSE.

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > Bempp::DiscreteBoundaryOperator< ValueType >::asDiscreteAcaBoundaryOperator ( double  eps = -1,
int  maximumRank = -1,
bool  interleave = false 
) const
virtual

Return a representation that can be cast to a DiscreteAcaBoundaryOperator.

The conversion only succeeds if all members of the DiscreteBoundaryOperator itself can be cast to a DiscreteAcaBoundaryOperator or if it is a linear combination of DiscreteOperators that can be cast to type DiscreteAcaBoundaryOperator. Operator compositions are not yet supported by this function.

Parameters
[in]epsThe $\epsilon$ parameter controlling the accuracy of H-matrix addition in case this operation is needed to construct the H-matrix representation of the operator. If eps is negative (default), it is set to the smaller of the values of $\epsilon$ specified during the creation of the two H-matrices to be added.
[in]maximumRankMaximum rank of blocks to be considered low rank when adding two H-Matrices. If maximumRank is negative, it is set to the larger of the maximum ranks specified during the creation of the two H-matrices to be added.
Returns
A pointer to a DiscreteBoundaryOperator object castable to DiscreteAcaBoundaryOperator.
Note
This function throws an exception if BEM++ has been compiled without AHMED.

Reimplemented in Bempp::DiscreteAcaBoundaryOperator< ValueType >, Bempp::DiscreteBlockedBoundaryOperator< ValueType >, Bempp::DiscreteSparseBoundaryOperator< ValueType >, Bempp::DiscreteBoundaryOperatorSum< ValueType >, and Bempp::ScaledDiscreteBoundaryOperator< ValueType >.

template<typename ValueType >
arma::Mat< ValueType > Bempp::DiscreteBoundaryOperator< ValueType >::asMatrix ( ) const
virtual
template<typename ValueType >
void Bempp::DiscreteBoundaryOperator< ValueType >::dump ( ) const
virtual

Write a textual representation of the operator to standard output.

The default implementation prints the matrix returned by asMatrix().

Reimplemented in Bempp::DiscreteSparseBoundaryOperator< ValueType >, Bempp::DiscreteDenseBoundaryOperator< ValueType >, and Bempp::DiscreteNullBoundaryOperator< ValueType >.

Friends And Related Function Documentation

template<typename RealType >
shared_ptr< DiscreteBoundaryOperator< std::complex< RealType > > > complexify ( const shared_ptr< const DiscreteBoundaryOperator< RealType > > &  op)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the "complexified" operator *op.

The returned operator has the same matrix representation as the original operator, but is able to act on vectors of complex numbers.

An exception is thrown if op is null.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > conjugate ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the conjugated operator *op.

An exception is thrown if op is null.

References Bempp::CONJUGATE.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > conjugateTranspose ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the conjugate-transposed operator *op.

An exception is thrown if op is null.

References Bempp::CONJUGATE_TRANSPOSE.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > mul ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op1,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op2 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing a composition of operators, (*op1) * (*op2).

An exception is thrown if op1 or op2 is null.

template<typename ValueType , typename ScalarType >
boost::enable_if< typename boost::mpl::has_key< boost::mpl::set< float, double, std::complex< float >, std::complex< double > >, ScalarType >, shared_ptr< DiscreteBoundaryOperator< ValueType > > >::type operator* ( ScalarType  scalar,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op multiplied by scalar.

An exception is thrown if op is null.

template<typename ValueType , typename ScalarType >
boost::enable_if< typename boost::mpl::has_key< boost::mpl::set< float, double, std::complex< float >, std::complex< double > >, ScalarType >, shared_ptr< DiscreteBoundaryOperator< ValueType > > >::type operator* ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op,
ScalarType  scalar 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op multiplied by scalar.

An exception is thrown if op is null.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > operator* ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op1,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op2 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing a composition of operators, (*op1) * (*op2).

An exception is thrown if op1 or op2 is null.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > operator+ ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op1,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op2 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the sum of the operands.

An exception is thrown if either op1 or op2 is null.

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > operator- ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the operand multiplied by -1.

An exception is thrown if op is null.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > operator- ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op1,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op2 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the difference of the operands.

An exception is thrown if either op1 or op2 is null.

template<typename ValueType , typename ScalarType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > operator/ ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op,
ScalarType  scalar 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the operator *op divided by scalar.

An exception is thrown if op is null.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > transpose ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the transposed operator *op.

An exception is thrown if op is null.

References Bempp::TRANSPOSE.

template<typename ValueType >
shared_ptr< DiscreteBoundaryOperator< ValueType > > transpose ( TranspositionMode  trans,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op 
)
related

Return a shared pointer to a DiscreteBoundaryOperator representing the transposed and/or conjugated (depending on the value of the argument trans) operator *op.

An exception is thrown if op is null.


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