BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Public Types | Public Member Functions | List of all members
Bempp::AssemblyOptions Class Reference

Options determining how weak-form assembly is done. More...

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

Public Types

enum  Value { AUTO = -1, NO = 0, YES = 1 }
 

Public Member Functions

void setMaxThreadCount (int maxThreadCount)
 Set the maximum number of threads used during the assembly. More...
 
BEMPP_DEPRECATED void switchToTbb (int maxThreadCount=AUTO)
 Set the maximum number of threads used during the assembly. More...
 
const ParallelizationOptionsparallelizationOptions () const
 Return current parallelization options.
 
void setVerbosityLevel (VerbosityLevel::Level level)
 Set the verbosity level. More...
 
VerbosityLevel::Level verbosityLevel () const
 Return the verbosity level.
 
void enableSingularIntegralCaching (bool value=true)
 Specify whether singular integrals are cached during weak-form assembly. More...
 
bool isSingularIntegralCachingEnabled () const
 Return whether singular integrals should be cached during weak-form assembly. More...
 
void enableSparseStorageOfLocalOperators (bool value=true)
 Specify whether discrete weak forms of local operators should be stored in sparse format. More...
 
bool isSparseStorageOfLocalOperatorsEnabled () const
 Return whether discrete weak forms of local operators should be stored in sparse format. More...
 
void BEMPP_DEPRECATED enableSparseStorageOfMassMatrices (bool value=true)
 Specify whether discrete weak forms of local operators should be stored in sparse format. More...
 
bool BEMPP_DEPRECATED isSparseStorageOfMassMatricesEnabled () const
 Return whether mass matrices should be stored in sparse format. More...
 
void enableJointAssembly (bool value=true)
 Enable or disable joint assembly of integral-operator superpositions. More...
 
bool isJointAssemblyEnabled () const
 Return whether joint assembly of integral-operator superpositions is enabled. More...
 
void enableBlasInQuadrature (Value value=AUTO)
 Specify whether BLAS matrix multiplication routines should be used during evaluation of elementary integrals. More...
 
Value isBlasEnabledInQuadrature () const
 Indicate whether BLAS matrix multiplication routines are used during evaluation of elementary integrals. More...
 
void makeQuadratureOrderUniformInEachCluster (bool value=true)
 Instruct the ACA assembler to use the same quadrature order for the evaluation of all regular integrals over pairs of test and trial elements belonging to a single block cluster. More...
 
bool isQuadratureOrderUniformInEachCluster () const
 Return whether the same quadrature order is used for all regular integrals over pairs of test and trial elements belonging to a single block cluster. More...
 

Assembly mode

enum  Mode { DENSE, ACA }
 Possible assembly modes for weak forms of boundary integral operators. More...
 
void switchToDenseMode ()
 Use dense-matrix representations of weak forms of boundary integral operators. More...
 
void switchToAcaMode (const AcaOptions &acaOptions)
 Use adaptive cross approximation (ACA) to obtain hierarchical-matrix representations of weak forms of boundary integral operators. More...
 
BEMPP_DEPRECATED void switchToDense ()
 Use dense-matrix representations of weak forms of boundary integral operators. More...
 
BEMPP_DEPRECATED void switchToAca (const AcaOptions &acaOptions)
 Use adaptive cross approximation (ACA) to obtain hierarchical-matrix representations of weak forms of boundary integral operators. More...
 
Mode assemblyMode () const
 Current assembly mode. More...
 
const AcaOptionsacaOptions () const
 Return the current adaptive cross approximation (ACA) settings. More...
 

Detailed Description

Options determining how weak-form assembly is done.

Member Enumeration Documentation

Possible assembly modes for weak forms of boundary integral operators.

Enumerator
DENSE 

Assemble dense matrices.

ACA 

Assemble hierarchical matrices using adaptive cross approximation (ACA).

Member Function Documentation

const AcaOptions & Bempp::AssemblyOptions::acaOptions ( ) const
AssemblyOptions::Mode Bempp::AssemblyOptions::assemblyMode ( ) const
void Bempp::AssemblyOptions::enableBlasInQuadrature ( Value  value = AUTO)

Specify whether BLAS matrix multiplication routines should be used during evaluation of elementary integrals.

If this option is set to AUTO (default), BLAS is used in the evaluation of integrals occurring in the weak forms of operators whose test or trial space is composed of quadratic or higher-order elements; for the remaining operators, BLAS routines are not used. You can force BLAS-based integration routines to be used always (or never) by setting this option to YES (or NO).

void Bempp::AssemblyOptions::enableJointAssembly ( bool  value = true)

Enable or disable joint assembly of integral-operator superpositions.

If value == true, the discrete weak forms of superpositions of elementary integral operators, such as 5. * SLP + 2. * DLP, will assembled together (for example, in a single run of ACA) and stored as a single H-matrix or dense matrix. By default joint assembly is disabled, which means that the discrete weak form of each elementary operator is stored separately.

void Bempp::AssemblyOptions::enableSingularIntegralCaching ( bool  value = true)

Specify whether singular integrals are cached during weak-form assembly.

If value == true, singular integrals are precalculated and stored in a cache before filling the matrix of the discretized weak form of a singular boundary operator. Otherwise these integrals are evaluated as needed during the assembly of the weak form.

By default, singular integral caching is enabled.

void Bempp::AssemblyOptions::enableSparseStorageOfLocalOperators ( bool  value = true)

Specify whether discrete weak forms of local operators should be stored in sparse format.

If value == true (default), discretized local operators are stored as sparse matrices, otherwise as dense matrices.

Referenced by enableSparseStorageOfMassMatrices().

void Bempp::AssemblyOptions::enableSparseStorageOfMassMatrices ( bool  value = true)

Specify whether discrete weak forms of local operators should be stored in sparse format.

Deprecated:
This function is deprecated and will be removed in a future version of BEM++. Use enableSparseStorageOfLocalOperators() instead.

References enableSparseStorageOfLocalOperators().

AssemblyOptions::Value Bempp::AssemblyOptions::isBlasEnabledInQuadrature ( ) const

Indicate whether BLAS matrix multiplication routines are used during evaluation of elementary integrals.

See enableBlasInQuadrature() for more information.

bool Bempp::AssemblyOptions::isJointAssemblyEnabled ( ) const

Return whether joint assembly of integral-operator superpositions is enabled.

See enableJointAssembly() for more information.

bool Bempp::AssemblyOptions::isQuadratureOrderUniformInEachCluster ( ) const

Return whether the same quadrature order is used for all regular integrals over pairs of test and trial elements belonging to a single block cluster.

See makeQuadratureOrderUniformInEachCluster() for more information.

bool Bempp::AssemblyOptions::isSingularIntegralCachingEnabled ( ) const

Return whether singular integrals should be cached during weak-form assembly.

See enableSingularIntegralCaching() for more information.

Referenced by Bempp::AbstractBoundaryOperator< BasisFunctionType_, ResultType_ >::collectOptionsDependentDataForAssemblerConstruction().

bool Bempp::AssemblyOptions::isSparseStorageOfLocalOperatorsEnabled ( ) const

Return whether discrete weak forms of local operators should be stored in sparse format.

See enableSparseStorageOfLocalOperators() for more information.

Referenced by isSparseStorageOfMassMatricesEnabled().

bool Bempp::AssemblyOptions::isSparseStorageOfMassMatricesEnabled ( ) const

Return whether mass matrices should be stored in sparse format.

See enableSparseStorageOfMassMatrices() for more information.

References isSparseStorageOfLocalOperatorsEnabled().

void Bempp::AssemblyOptions::makeQuadratureOrderUniformInEachCluster ( bool  value = true)

Instruct the ACA assembler to use the same quadrature order for the evaluation of all regular integrals over pairs of test and trial elements belonging to a single block cluster.

Use of different quadrature rules to evaluate regular integrals making up the entries of a single block of an H-matrix may lead to an increase in the block rank, and hence memory consumption, if these quadrature rules are not accurate enough (if their error is greater than the ACA tolerance $\epsilon$). However, evaluation of the entries of large blocks using a uniform quadrature rule may be unnecessarily time-consuming.

If this option is set to false and one uses a quadrature strategy that takes the interlement distance into account when selecting quadrature order, the true interelement distance is calculated and used for each pair of elements. Otherwise, instead of the distance between elements the quadrature-rule selector is passed the distance between the clusters to which these elements belong.

By default, this option is set to true.

void Bempp::AssemblyOptions::setMaxThreadCount ( int  maxThreadCount)

Set the maximum number of threads used during the assembly.

maxThreadCount must be a positive number or AUTO. In the latter case the number of threads is determined automatically.

Referenced by switchToTbb().

void Bempp::AssemblyOptions::setVerbosityLevel ( VerbosityLevel::Level  level)

Set the verbosity level.

This setting determines the amount of information printed out by functions from BEM++.

Referenced by Bempp::SyntheticIntegralOperator< BasisFunctionType_, ResultType_ >::getContextsForInternalAndAuxiliaryOperators().

void Bempp::AssemblyOptions::switchToAca ( const AcaOptions acaOptions)

Use adaptive cross approximation (ACA) to obtain hierarchical-matrix representations of weak forms of boundary integral operators.

Deprecated:
Use switchToAcaMode() instead.

References switchToAcaMode().

void Bempp::AssemblyOptions::switchToAcaMode ( const AcaOptions acaOptions)

Use adaptive cross approximation (ACA) to obtain hierarchical-matrix representations of weak forms of boundary integral operators.

Parameters
[in]acaOptionsParameters influencing the ACA algorithm.

References ACA, acaOptions(), Bempp::AcaOptions::globalAssemblyBeforeCompression, Bempp::AcaOptions::mode, and Bempp::AcaOptions::reactionToUnsupportedMode.

Referenced by Bempp::SyntheticIntegralOperator< BasisFunctionType_, ResultType_ >::getContextsForInternalAndAuxiliaryOperators(), and switchToAca().

void Bempp::AssemblyOptions::switchToDense ( )

Use dense-matrix representations of weak forms of boundary integral operators.

Deprecated:
Use switchToDenseMode() instead.

References switchToDenseMode().

void Bempp::AssemblyOptions::switchToDenseMode ( )

Use dense-matrix representations of weak forms of boundary integral operators.

This is the default assembly mode.

References DENSE.

Referenced by switchToDense().

void Bempp::AssemblyOptions::switchToTbb ( int  maxThreadCount = AUTO)

Set the maximum number of threads used during the assembly.

Deprecated:
Use setMaxThreadCount() instead.

References setMaxThreadCount().


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