OpenBEM
Open-source framework for electromagnetic simulation with the boundary element method.
Loading...
Searching...
No Matches
tmfie.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_TMFIE_H
19#define BEM_RWG_INTEQ_TMFIE_H
20
21#include "types.hpp"
22#include "constants.hpp"
23#include "materials.hpp"
24
26
27#include "matrix/base.hpp"
29
31
35
38
40
44
45
46namespace bem::rwg
47{
48
58template <typename MatrixType = EigenDenseMatrix<Complex>>
59class Tmfie: public IntegralEquationBase<MatrixType>
60{
61
63 using base::base;
64
65public:
66
75 )
76 {
77 op_T_ = op_T;
78 op_K_ = op_K;
79 return;
80 };
81
82
89 MatrixType j_matrix(
90 const Float f,
91 const Material& material
92 )
93 {
94 Complex k = material.k(f);
95
96 MatrixType K;
97 assembler_.assemble(K, op_K_, k);
98
99 K.add_ax(id_matrix(), half);
100
101 return K;
102 };
103
104
111 MatrixType m_matrix(
112 const Float f,
113 const Material& material
114 )
115 {
116 Complex k = material.k(f);
117 Complex eps = material.eps();
118 Complex sigma = material.sigma();
119
120 MatrixType T;
121 assembler_.assemble(T, op_T_, k);
122
123 T.scale(-(J * two_pi * f * eps + sigma));
124
125 return T;
126 };
127
128
133 MatrixType id_matrix()
134 {
136 MatrixType Ir;
137 assembler_.assemble(Ir, op_Ir, 0);
138 if (base::flip_normals_)
139 Ir.scale(-one);
140 return Ir;
141 };
142
143
151 MatrixType exc_matrix(
152 const Float f,
153 const Material& material,
155 )
156 {
157 std::vector<Index> obs_elems_vec (base::elem_pairs_.cols());
158 EigRowVec<Index>::Map(&obs_elems_vec[0], base::elem_pairs_.row(0).size()) = base::elem_pairs_.row(0);
159
160 std::sort(obs_elems_vec.begin(), obs_elems_vec.end());
161 obs_elems_vec.erase(std::unique(obs_elems_vec.begin(), obs_elems_vec.end()), obs_elems_vec.end());
162
165
166 Complex k = material.k(f);
167 MatrixType inc (base::obs_mesh_.num_edges(), exc.num_excitations());
168 exc_assembler.assemble(inc, exc, k);
169
170 return inc;
171 };
172
173
181 MatrixType j_projector(
182 const Float f,
183 const Material& material,
184 const PointCloud<3>& points
185 )
186 {
187 Complex k = material.k(f);
188
189 EdgeProjectorAssembler<3> proj_assembler (points, base::src_mesh_);
190
191 MatrixType K;
192 proj_assembler.assemble(K, proj_K_, k);
193
194 return K;
195 };
196
197
205 MatrixType m_projector(
206 const Float f,
207 const Material& material,
208 const PointCloud<3>& points
209 )
210 {
211 Complex k = material.k(f);
212 Complex eps = material.eps();
213 Complex sigma = material.sigma();
214
215 EdgeProjectorAssembler<3> proj_assembler (points, base::src_mesh_);
216
217 MatrixType T;
218 proj_assembler.assemble(T, proj_T_, k);
219
220 T.scale(-(J * two_pi * f * eps + sigma));
221
222 return T;
223 };
224
225
226protected:
227
230
233
234 EdgeOperatorAssembler assembler_ = EdgeOperatorAssembler(base::obs_mesh_, base::src_mesh_, base::elem_pairs_);
235
236};
237
242}
243
244#endif
Class defining a general material with a constant (zero or non-zero) electrical conductivity and real...
Definition materials.hpp:43
virtual Complex eps() const
Returns the permittivity.
Definition materials.hpp:64
virtual Complex k(const Float f) const
Returns the frequency-dependent complex wavenumber.
virtual Float sigma() const
Returns the electrical conductivity.
Class for generating excitation matrices for RWG-based BEM systems with RWG testing functions.
Class for generating operator matrices for RWG observation and source functions.
Base class defining an RWG-based integral equation.
Definition base.hpp:44
void assemble(MatrixBase< Complex > &mat, OperatorType op, const Complex k)
Assembles the operator matrix for a given operator object and source and observation meshes.
Definition base.hpp:115
Class defining the RWG-based TMFIE.
Definition tmfie.hpp:60
MatrixType m_matrix(const Float f, const Material &material)
Returns the operator matrix associated with the magnetic surface current density.
Definition tmfie.hpp:111
MatrixType j_matrix(const Float f, const Material &material)
Returns the operator matrix associated with the electric surface current density.
Definition tmfie.hpp:89
MatrixType m_projector(const Float f, const Material &material, const PointCloud< 3 > &points)
Returns the projector matrix associated with the magnetic surface current density.
Definition tmfie.hpp:205
void set_operators(const VectorHypersingularOp<> &op_T, const VectorDoubleLayerPvOp<> &op_K)
Sets custom operators for the Tmfie object.
Definition tmfie.hpp:72
MatrixType exc_matrix(const Float f, const Material &material, ExcitationBase< 3 > &exc)
Returns the excitation matrix for a given excitation operator.
Definition tmfie.hpp:151
MatrixType id_matrix()
Returns the identity operator matrix associated with the electric surface current density.
Definition tmfie.hpp:133
MatrixType j_projector(const Float f, const Material &material, const PointCloud< 3 > &points)
Returns the projector matrix associated with the electric surface current density.
Definition tmfie.hpp:181
const Complex J
Imaginary unit.
Definition constants.hpp:40
Eigen::Matrix< T, 1, Eigen::Dynamic > EigRowVec
Dynamic-size row vector containing type T.
Definition types.hpp:90
double Float
Floating point number.
Definition types.hpp:47
std::complex< Float > Complex
Complex floating point number.
Definition types.hpp:51
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.