OpenBEM
Open-source framework for electromagnetic simulation with the boundary element method.
Loading...
Searching...
No Matches
base.hpp
Go to the documentation of this file.
1// OpenBEM - Copyright (C) 2026 Shashwat Sharma
2
3// This file is part of OpenBEM.
4
5// OpenBEM is free software: you can redistribute it and/or modify it under the terms of the
6// GNU General Public License as published by the Free Software Foundation, either version 3
7// of the License, or (at your option) any later version.
8
9// You should have received a copy of the GNU General Public License along with OpenBEM.
10// If not, see <https://www.gnu.org/licenses/>.
11
12
18#ifndef BEM_RWG_INTEQ_BASE_H
19#define BEM_RWG_INTEQ_BASE_H
20
21#include "types.hpp"
23
24#include "matrix/base.hpp"
26
28
29
30namespace bem::rwg
31{
32
42template <typename MatrixType = EigenDenseMatrix<Complex>>
44{
45
46static_assert(
47 std::is_base_of<MatrixBase<Complex>, MatrixType>::value,
48 "IntegralEquationBase: MatrixType must be derived from MatrixBase<Complex>."
49 );
50
51public:
52
65 const bool flip_normals = false
66 ):
67 obs_mesh_(obs_mesh),
68 src_mesh_(src_mesh),
69 elem_pairs_(elem_pairs),
70 flip_normals_(flip_normals) {};
71
72
82 const bool flip_normals = false
83 ):
84 obs_mesh_(obs_mesh),
85 src_mesh_(src_mesh),
86 elem_pairs_(OperatorAssemblerBase<3, 3>::make_pairs(obs_mesh, src_mesh)),
87 flip_normals_(flip_normals) {};
88
89
98 const TriangleMesh<3>& mesh,
100 const bool flip_normals = false
101 ):
102 obs_mesh_(mesh),
103 src_mesh_(mesh),
104 elem_pairs_(elem_pairs),
105 flip_normals_(flip_normals) {};
106
107
114 const TriangleMesh<3>& mesh,
115 const bool flip_normals = false
116 ):
117 obs_mesh_(mesh),
118 src_mesh_(mesh),
119 elem_pairs_(OperatorAssemblerBase<3, 3>::make_pairs(mesh, mesh)),
120 flip_normals_(flip_normals) {};
121
122
123protected:
124
125 const TriangleMesh<3>& obs_mesh_;
126 const TriangleMesh<3>& src_mesh_;
127 const EigMatNX<Index, 2> elem_pairs_;
128 const bool flip_normals_ = false;
129
130};
131
136}
137
138#endif
Base class defining an RWG-based integral equation.
Definition base.hpp:44
IntegralEquationBase(const TriangleMesh< 3 > &obs_mesh, const TriangleMesh< 3 > &src_mesh, ConstEigRef< EigMatNX< Index, 2 > > elem_pairs, const bool flip_normals=false)
Constructs an IntegralEquationBase object.
Definition base.hpp:61
IntegralEquationBase(const TriangleMesh< 3 > &mesh, const bool flip_normals=false)
Constructs an IntegralEquationBase object.
Definition base.hpp:113
IntegralEquationBase(const TriangleMesh< 3 > &mesh, ConstEigRef< EigMatNX< Index, 2 > > elem_pairs, const bool flip_normals=false)
Constructs an IntegralEquationBase object.
Definition base.hpp:97
IntegralEquationBase(const TriangleMesh< 3 > &obs_mesh, const TriangleMesh< 3 > &src_mesh, const bool flip_normals=false)
Constructs an IntegralEquationBase object.
Definition base.hpp:79
Base class for generating RWG-based BEM operator matrices.
Definition base.hpp:46
const Eigen::Ref< const EigObj > ConstEigRef
Read-only reference to an Eigen object.
Definition types.hpp:98
Eigen::Matrix< T, N, 1 > EigColVecN
Fixed-size column vector of size N containing type T.
Definition types.hpp:86
Namespace for RWG-based BEM functionality.