OpenBEM
Open-source framework for electromagnetic simulation with the boundary element method.
Loading...
Searching...
No Matches
bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof > Class Template Referenceabstract

Base class for generating RWG-based BEM operator matrices. More...

#include <base.hpp>

Detailed Description

template<uint8_t obs_num_dof, uint8_t src_num_dof>
class bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >

Base class for generating RWG-based BEM operator matrices.

Template Parameters
obs_num_dof- Number of degrees of freedom associated with each observation triangle.
src_num_dof- Number of degrees of freedom associated with each source triangle.

Definition at line 45 of file base.hpp.

Constructor & Destructor Documentation

◆ OperatorAssemblerBase() [1/4]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::OperatorAssemblerBase ( const TriangleMesh< 3 > &  obs_mesh,
const TriangleMesh< 3 > &  src_mesh,
ConstEigRef< EigMatNX< Index, 2 > >  elem_pairs 
)
inline

Constructs an OperatorAssemblerBase for given observation and source meshes.

Parameters
[in]obs_mesh- Observation triangle mesh for which the operator matrix is to be assembled.
[in]src_mesh- Source triangle mesh for which the operator matrix is to be assembled.
[in]elem_pairs- Observation (first row) and source (second row) triangle index pairs for which the operator matrix is to be assembled.

Definition at line 60 of file base.hpp.

◆ OperatorAssemblerBase() [2/4]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::OperatorAssemblerBase ( const TriangleMesh< 3 > &  obs_mesh,
const TriangleMesh< 3 > &  src_mesh 
)
inline

Constructs an OperatorAssemblerBase for given observation and source meshes.

Parameters
[in]obs_mesh- Observation triangle mesh for which the operator matrix is to be assembled.
[in]src_mesh- Source triangle mesh for which the operator matrix is to be assembled.

Definition at line 75 of file base.hpp.

◆ OperatorAssemblerBase() [3/4]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::OperatorAssemblerBase ( const TriangleMesh< 3 > &  mesh,
ConstEigRef< EigMatNX< Index, 2 > >  elem_pairs 
)
inline

Constructs an OperatorAssemblerBase for a given mesh.

Parameters
[in]mesh- Triangle mesh for which the operator matrix is to be assembled.
[in]elem_pairs- Observation (first row) and source (second row) triangle index pairs for which the operator matrix is to be assembled.

Definition at line 90 of file base.hpp.

◆ OperatorAssemblerBase() [4/4]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::OperatorAssemblerBase ( const TriangleMesh< 3 > &  mesh)
inline

Constructs an OperatorAssemblerBase for a given mesh.

Parameters
[in]mesh- Triangle mesh for which the operator matrix is to be assembled.

Definition at line 100 of file base.hpp.

Member Function Documentation

◆ assemble()

template<uint8_t obs_num_dof, uint8_t src_num_dof>
template<typename OperatorType >
void bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::assemble ( MatrixBase< Complex > &  mat,
OperatorType  op,
const Complex  k 
)
inline

Assembles the operator matrix for a given operator object and source and observation meshes.

Parameters
[out]mat- Matrix to store the assembled operator coefficients, with columns corresponding to source degrees of freedom, and rows corresponding to observation degrees of freedom.
[in]op- Operator object that computes the coefficients to be assembled into mat; must inherit from OperatorBase<obs_num_dof, src_num_dof>.
[in]k- Complex wavenumber.

Definition at line 115 of file base.hpp.

◆ prep_matrix()

template<uint8_t obs_num_dof, uint8_t src_num_dof>
virtual void bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::prep_matrix ( MatrixBase< Complex > &  mat)
pure virtual

Prepares the matrix for assembly (e.g., resizing and preallocation).

Parameters
[out]mat- Matrix to store the assembled operator coefficients, with columns corresponding to source degrees of freedom, and rows corresponding to observation degrees of freedom.

Implemented in bem::rwg::EdgeOperatorAssembler, bem::rwg::FaceOperatorAssembler, bem::rwg::FaceEdgeOperatorAssembler, bem::rwg::EdgeFaceOperatorAssembler, and bem::rwg::VectorOperatorsAssembler.

◆ fill_matrix()

template<uint8_t obs_num_dof, uint8_t src_num_dof>
virtual void bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::fill_matrix ( MatrixBase< Complex > &  mat,
ConstEigRef< EigColVecN< Index, 2 > >  elem_pair,
ConstEigRef< EigMatMN< Complex, obs_num_dof, src_num_dof > >  values 
)
pure virtual

Fills operator values in the matrix based on source and observation meshes and degrees of freedom.

Parameters
[out]mat- Matrix to store the assembled operator coefficients, with columns corresponding to source degrees of freedom, and rows corresponding to observation degrees of freedom.
[in]elem_pair- Observation (first entry) and source (second entry) triangle index pair.
[in]values- Operator values for each pair of observation and source degrees of freedom.

◆ make_pairs() [1/2]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
static EigMatNX< Index, 2 > bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::make_pairs ( const TriangleMesh< 3 > &  obs_mesh,
const TriangleMesh< 3 > &  src_mesh 
)
inlinestatic

Generates all possible pairs of triangle indices for given observation and source triangle meshes.

Parameters
[in]obs_mesh- Observation triangle mesh.
[in]src_mesh- Source triangle mesh.
Returns
All possible triangle index pairs, with observation indices in the first row, and source indices in the second row.

Definition at line 175 of file base.hpp.

◆ make_pairs() [2/2]

template<uint8_t obs_num_dof, uint8_t src_num_dof>
static EigMatNX< Index, 2 > bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::make_pairs ( ConstEigRef< EigRowVec< Index > >  obs_elems,
ConstEigRef< EigRowVec< Index > >  src_elems 
)
inlinestatic

Makes all possible pairs of triangle indices for given observation and source triangle indices.

Parameters
[in]obs_elems- Observation triangle indices.
[in]src_elems- Source triangle indices.
Returns
All possible triangle index pairs, with observation indices in the first row, and source indices in the second row.

Definition at line 200 of file base.hpp.

Member Data Documentation

◆ obs_mesh_

template<uint8_t obs_num_dof, uint8_t src_num_dof>
const TriangleMesh<3>& bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::obs_mesh_
protected

Definition at line 223 of file base.hpp.

◆ src_mesh_

template<uint8_t obs_num_dof, uint8_t src_num_dof>
const TriangleMesh<3>& bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::src_mesh_
protected

Definition at line 224 of file base.hpp.

◆ elem_pairs_

template<uint8_t obs_num_dof, uint8_t src_num_dof>
const EigMatNX<Index, 2> bem::rwg::OperatorAssemblerBase< obs_num_dof, src_num_dof >::elem_pairs_
protected

Definition at line 225 of file base.hpp.


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