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

Discrete linear operator stored as a H-matrix. More...

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

Inheritance diagram for Bempp::DiscreteAcaBoundaryOperator< ValueType >:
Bempp::DiscreteBoundaryOperator< ValueType >

Public Types

typedef Fiber::ScalarTraits
< ValueType >::RealType 
CoordinateType
 
typedef AhmedDofWrapper
< CoordinateType > 
AhmedDofType
 
typedef bbxbemblcluster
< AhmedDofType, AhmedDofType
AhmedBemBlcluster
 
typedef mblock< typename
AhmedTypeTraits< ValueType >
::Type > 
AhmedMblock
 
typedef boost::shared_array
< AhmedMblock * > 
AhmedMblockArray
 
typedef boost::shared_array
< AhmedMblock * > 
AhmedConstMblockArray
 

Public Member Functions

 DiscreteAcaBoundaryOperator (unsigned int rowCount, unsigned int columnCount, double epsUsedInAssembly, int maximumRankUsedInAssembly, int symmetry, const shared_ptr< const AhmedBemBlcluster > &blockCluster_, const AhmedMblockArray &blocks_, const IndexPermutation &domainPermutation_, const IndexPermutation &rangePermutation_, const ParallelizationOptions &parallelizationOptions_, const std::vector< AhmedConstMblockArray > &sharedBlocks_=std::vector< AhmedConstMblockArray >())
 Constructor. More...
 
 DiscreteAcaBoundaryOperator (unsigned int rowCount, unsigned int columnCount, int maximumRankUsedInAssembly, int symmetry, std::auto_ptr< const AhmedBemBlcluster > blockCluster_, AhmedMblockArray blocks_, const IndexPermutation &domainPermutation_, const IndexPermutation &rangePermutation_, const ParallelizationOptions &parallelizationOptions_, const std::vector< AhmedConstMblockArray > &sharedBlocks_=std::vector< AhmedConstMblockArray >()) BEMPP_DEPRECATED
 Constructor. More...
 
BEMPP_DEPRECATED DiscreteAcaBoundaryOperator (unsigned int rowCount, unsigned int columnCount, int maximumRankUsedInAssembly, int symmetry, const shared_ptr< const AhmedBemBlcluster > &blockCluster_, const AhmedMblockArray &blocks_, const IndexPermutation &domainPermutation_, const IndexPermutation &rangePermutation_, const ParallelizationOptions &parallelizationOptions_, const std::vector< AhmedConstMblockArray > &sharedBlocks_=std::vector< AhmedConstMblockArray >())
 Constructor. More...
 
virtual arma::Mat< ValueType > asMatrix () const
 Matrix representation of the operator. More...
 
virtual unsigned int rowCount () const
 Number of rows of the operator.
 
virtual unsigned int columnCount () const
 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
 Add a subblock of this operator to a matrix. More...
 
virtual shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
asDiscreteAcaBoundaryOperator (double eps=-1, int maximumRank=-1, bool interleave=false) const
 Return a representation that can be cast to a DiscreteAcaBoundaryOperator. More...
 
void makeAllMblocksDense ()
 Uncompress all blocks of the H-matrix and store them as dense matrices. More...
 
int maximumRank () const
 Return the upper bound for the rank of low-rank mblocks specified during H-matrix construction.
 
double eps () const
 Return the value of the epsilon parameter specified during H-matrix construction.
 
int actualMaximumRank () const
 Return the actual maximum rank of low-rank mblocks.
 
int symmetry () const
 Return a flag describing the symmetry of this operator.
 
shared_ptr< const
AhmedBemBlcluster > 
blockCluster () const
 Return the block cluster.
 
AhmedMblockArray blocks () const
 Return the mblocks making up the H-matrix represented by this operator. More...
 
size_t blockCount () const
 Return the number of mblocks making up this operator.
 
const IndexPermutationdomainPermutation () const
 Return the domain index permutation.
 
const IndexPermutationrangePermutation () const
 Return the range index permutation.
 
const ParallelizationOptionsparallelizationOptions () const
 Return the parallelization options used in the matrix-vector multiply.
 
std::vector
< AhmedConstMblockArray > 
sharedBlocks () const
 Return the vector of mblock arrays that this operator implicitly depends on.
 
virtual Teuchos::RCP< const
Thyra::VectorSpaceBase
< ValueType > > 
domain () const
 
virtual Teuchos::RCP< const
Thyra::VectorSpaceBase
< ValueType > > 
range () const
 
- Public Member Functions inherited from Bempp::DiscreteBoundaryOperator< ValueType >
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 void dump () const
 Write a textual representation of the operator to standard output. More...
 

Static Public Member Functions

static const
DiscreteAcaBoundaryOperator
< ValueType > & 
castToAca (const DiscreteBoundaryOperator< ValueType > &discreteOperator)
 Downcast a reference to a DiscreteBoundaryOperator object to DiscreteAcaBoundaryOperator. More...
 
static shared_ptr< const
DiscreteAcaBoundaryOperator
< ValueType > > 
castToAca (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &discreteOperator)
 Downcast a shared pointer to a DiscreteBoundaryOperator object to a shared pointer to a DiscreteAcaBoundaryOperator. More...
 

Protected Member Functions

virtual bool opSupportedImpl (Thyra::EOpTransp M_trans) const
 
- Protected Member Functions inherited from Bempp::DiscreteBoundaryOperator< ValueType >
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
 

Friends

class AcaApproximateLuInverse< ValueType >
 
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
acaOperatorSum (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2, double eps, int maximumRank)
 
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
scaledAcaOperator (const ValueType &multiplier, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
scaledAcaOperator (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op, const ValueType &multiplier)
 

Related Functions

(Note that these are not member functions.)

template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
acaOperatorSum (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op1, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op2, double eps, int maximumRank)
 Add two discrete boundary operators stored as H-matrices. More...
 
template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
scaledAcaOperator (const ValueType &multiplier, const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op)
 Multiply the H-matrix representation of a discrete boundary operator by a scalar and wrap the result in a new discrete boundary operator. More...
 
template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
scaledAcaOperator (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op, const ValueType &multiplier)
 Multiply the H-matrix representation of a discrete boundary operator by a scalar and wrap the result in a new discrete boundary operator. More...
 
template<typename ValueType >
shared_ptr< const
DiscreteBoundaryOperator
< ValueType > > 
acaOperatorApproximateLuInverse (const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &op, double delta)
 LU inverse of a discrete boundary operator stored as a H-matrix. More...
 

Detailed Description

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

Discrete linear operator stored as a H-matrix.

Constructor & Destructor Documentation

template<typename ValueType >
Bempp::DiscreteAcaBoundaryOperator< ValueType >::DiscreteAcaBoundaryOperator ( unsigned int  rowCount,
unsigned int  columnCount,
double  epsUsedInAssembly,
int  maximumRankUsedInAssembly,
int  symmetry,
const shared_ptr< const AhmedBemBlcluster > &  blockCluster_,
const AhmedMblockArray &  blocks_,
const IndexPermutation domainPermutation_,
const IndexPermutation rangePermutation_,
const ParallelizationOptions parallelizationOptions_,
const std::vector< AhmedConstMblockArray > &  sharedBlocks_ = std::vector<AhmedConstMblockArray>() 
)

Constructor.

Parameters
[in]rowCountNumber of rows.
[in]columnCountNumber of columns.
[in]epsUsedInAssemblyThe epsilon parameter used during assembly of the H-matrix.
[in]maximumRankUsedInAssemblyThe limit on block rank used during assembly of the H-matrix. This is used in particular when generating the approximate H-matrix LU.
[in]symmetryH-matrix symmetry. Can be any combination of the flags defined in the Symmetry enumeration type.
[in]blockCluster_Block cluster defining the structure of the H-matrix.
[in]blocks_Array containing the H-matrix blocks.
[in]domainPermutation_Mapping from original to permuted column indices.
[in]rangePermutation_Mapping from original to permuted row indices.
[in]parallelizationOptions_Options determining the maximum number of threads used in the apply() routine for the H-matrix-vector product.
[in]sharedBlocks_Vector of arrays of mblocks on which this operator implicitly depends and which therefore must stay alive for the lifetime of this operator. Useful for constructing ACA operators that combine mblocks of several other operators.
Note
Currently the apply() routine is only parallelized for non-Hermitian H-matrices.
template<typename ValueType >
Bempp::DiscreteAcaBoundaryOperator< ValueType >::DiscreteAcaBoundaryOperator ( unsigned int  rowCount,
unsigned int  columnCount,
int  maximumRankUsedInAssembly,
int  symmetry,
std::auto_ptr< const AhmedBemBlcluster >  blockCluster_,
AhmedMblockArray  blocks_,
const IndexPermutation domainPermutation_,
const IndexPermutation rangePermutation_,
const ParallelizationOptions parallelizationOptions_,
const std::vector< AhmedConstMblockArray > &  sharedBlocks_ = std::vector<AhmedConstMblockArray>() 
)

Constructor.

Parameters
[in]rowCountNumber of rows.
[in]columnCountNumber of columns.
[in]maximumRankUsedInAssemblyThe limit on block rank used during assembly of the H-matrix. This is used in particular when generating the approximate H-matrix LU.
[in]symmetryH-matrix symmetry. Can be any combination of the flags defined in the Symmetry enumeration type.
[in]blockCluster_Block cluster defining the structure of the H-matrix.
[in]blocks_Array containing the H-matrix blocks.
[in]domainPermutation_Mapping from original to permuted column indices.
[in]rangePermutation_Mapping from original to permuted row indices.
[in]parallelizationOptions_Options determining the maximum number of threads used in the apply() routine for the H-matrix-vector product.
[in]sharedBlocks_Vector of arrays of mblocks on which this operator implicitly depends and which therefore must stay alive for the lifetime of this operator. Useful for constructing ACA operators that combine mblocks of several other operators.
Note
Currently the apply() routine is only parallelized for non-Hermitian H-matrices.
Deprecated:
This constructor is deprecated. Use the non-deprecated constructor.
template<typename ValueType >
Bempp::DiscreteAcaBoundaryOperator< ValueType >::DiscreteAcaBoundaryOperator ( unsigned int  rowCount,
unsigned int  columnCount,
int  maximumRankUsedInAssembly,
int  symmetry,
const shared_ptr< const AhmedBemBlcluster > &  blockCluster_,
const AhmedMblockArray &  blocks_,
const IndexPermutation domainPermutation_,
const IndexPermutation rangePermutation_,
const ParallelizationOptions parallelizationOptions_,
const std::vector< AhmedConstMblockArray > &  sharedBlocks_ = std::vector<AhmedConstMblockArray>() 
)

Constructor.

Parameters
[in]rowCountNumber of rows.
[in]columnCountNumber of columns.
[in]maximumRankUsedInAssemblyThe limit on block rank used during assembly of the H-matrix. This is used in particular when generating the approximate H-matrix LU.
[in]symmetryH-matrix symmetry. Can be any combination of the flags defined in the Symmetry enumeration type.
[in]blockCluster_Block cluster defining the structure of the H-matrix.
[in]blocks_Array containing the H-matrix blocks.
[in]domainPermutation_Mapping from original to permuted column indices.
[in]rangePermutation_Mapping from original to permuted row indices.
[in]parallelizationOptions_Options determining the maximum number of threads used in the apply() routine for the H-matrix-vector product.
[in]sharedBlocks_Vector of arrays of mblocks on which this operator implicitly depends and which therefore must stay alive for the lifetime of this operator. Useful for constructing ACA operators that combine mblocks of several other operators.
Note
Currently the apply() routine is only parallelized for non-Hermitian H-matrices.
Deprecated:
This constructor is deprecated. Use the non-deprecated constructor.

Member Function Documentation

template<typename ValueType>
void Bempp::DiscreteAcaBoundaryOperator< ValueType >::addBlock ( const std::vector< int > &  rows,
const std::vector< int > &  cols,
const ValueType  alpha,
arma::Mat< ValueType > &  block 
) const
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.

Implements Bempp::DiscreteBoundaryOperator< ValueType >.

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > Bempp::DiscreteAcaBoundaryOperator< 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 from Bempp::DiscreteBoundaryOperator< ValueType >.

template<typename ValueType >
arma::Mat< ValueType > Bempp::DiscreteAcaBoundaryOperator< ValueType >::asMatrix ( ) const
virtual

Matrix representation of the operator.

The default implementation is slow and should be overridden where possible.

Reimplemented from Bempp::DiscreteBoundaryOperator< ValueType >.

template<typename ValueType >
DiscreteAcaBoundaryOperator< ValueType >::AhmedMblockArray Bempp::DiscreteAcaBoundaryOperator< ValueType >::blocks ( ) const

Return the mblocks making up the H-matrix represented by this operator.

Note
This function returns a shared array of pointers to non-constant blocks. However, you must not modify it! This is just a workaround for AHMED's lack of const-correctness.
template<typename ValueType>
const DiscreteAcaBoundaryOperator< ValueType > & Bempp::DiscreteAcaBoundaryOperator< ValueType >::castToAca ( const DiscreteBoundaryOperator< ValueType > &  discreteOperator)
static
template<typename ValueType>
shared_ptr< const DiscreteAcaBoundaryOperator< ValueType > > Bempp::DiscreteAcaBoundaryOperator< ValueType >::castToAca ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  discreteOperator)
static

Downcast a shared pointer to a DiscreteBoundaryOperator object to a shared pointer to a DiscreteAcaBoundaryOperator.

If the object referenced by discreteOperator is not in fact a DiscreteAcaBoundaryOperator, a std::bad_cast exception is thrown.

template<typename ValueType >
void Bempp::DiscreteAcaBoundaryOperator< ValueType >::makeAllMblocksDense ( )

Uncompress all blocks of the H-matrix and store them as dense matrices.

Sometimes useful for debugging.

Friends And Related Function Documentation

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > acaOperatorApproximateLuInverse ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op,
double  delta 
)
related

LU inverse of a discrete boundary operator stored as a H-matrix.

Parameters
[in]opDiscrete boundary operator for which to compute the LU inverse.
[in]deltaApproximation accuracy of the inverse.
Returns
A shared pointer to a newly allocated discrete boundary operator representing the (approximate) LU inverse of op and stored as an (approximate) LU decomposition of op.

References Bempp::DiscreteAcaBoundaryOperator< ValueType >::castToAca().

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > acaOperatorSum ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op1,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op2,
double  eps,
int  maximumRank 
)
related

Add two discrete boundary operators stored as H-matrices.

A std::bad_cast exception is thrown if the input operators can not be cast to DiscreteAcaBoundaryOperator.

Parameters
[in]op1First operand.
[in]op2Second operand.
[in]eps???
Todo:
look into M. Bebendorf's book.
Parameters
[in]maximumRankMaximum rank of blocks that should be considered low-rank in the H-matrix to be constructed.
Returns
A shared pointer to a newly allocated discrete boundary operator representing the sum of the operands op1 and op2 stored as a single H-matrix.

References Bempp::DiscreteAcaBoundaryOperator< ValueType >::castToAca().

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > scaledAcaOperator ( const ValueType &  multiplier,
const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op 
)
related

Multiply the H-matrix representation of a discrete boundary operator by a scalar and wrap the result in a new discrete boundary operator.

A std::bad_cast exception is thrown if the input operator can not be cast to DiscreteAcaBoundaryOperator

Parameters
[in]multiplierScalar multiplier.
[in]opDiscrete boundary operator to be multiplied.
Returns
A shared pointer to a newly allocated discrete boundary operator representing the operand op multiplied by multiplier and stored as a H-matrix.

References Bempp::DiscreteAcaBoundaryOperator< ValueType >::castToAca().

template<typename ValueType >
shared_ptr< const DiscreteBoundaryOperator< ValueType > > scaledAcaOperator ( const shared_ptr< const DiscreteBoundaryOperator< ValueType > > &  op,
const ValueType &  multiplier 
)
related

Multiply the H-matrix representation of a discrete boundary operator by a scalar and wrap the result in a new discrete boundary operator.

A std::bad_cast exception is thrown if the input operator can not be cast to DiscreteAcaBoundaryOperator

Parameters
[in]opDiscrete boundary operator to be multiplied.
[in]multiplierScalar multiplier.
Returns
A shared pointer to a newly allocated discrete boundary operator representing the operand op multiplied by multiplier and stored as a H-matrix.

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