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_OPS_BASE_H
19#define BEM_RWG_OPS_BASE_H
20
21#include "types.hpp"
24
25
26namespace bem::rwg
27{
28
39template <uint8_t obs_num_dof, uint8_t src_num_dof>
41{
42
43 static_assert((obs_num_dof > 0), "OperatorBase: `obs_num_dof` must be greater than 0.");
44 static_assert((src_num_dof > 0), "OperatorBase: `src_num_dof` must be greater than 0.");
45
46public:
47
59 const Complex k,
60 const Triangle<3>& obs_tri,
61 const Triangle<3>& src_tri
62 ) = 0;
63
64
76 const Triangle<3>& obs_tri,
77 const Triangle<3>& src_tri
78 )
79 {
80 EigMatMN<Float, 3, 3> local_uvw = src_tri.local_coordinate_basis();
81 EigColVecN<Float, 3> local_origin = src_tri.local_origin();
82
83 obs_tri_local.set_v(
85 obs_tri.v(), local_origin, local_uvw
86 ),
87 obs_tri.edge_polarities()
88 );
89 src_tri_local.set_v(
91 src_tri.v(), local_origin, local_uvw
92 ).topRows(2),
93 src_tri.edge_polarities()
94 );
95
96 return;
97 };
98
99};
100
105}
106
107#endif
Geometry operations class.
Base class for RWG-based BEM operators.
Definition base.hpp:41
static void transform_coordinates(Triangle< 3 > &obs_tri_local, Triangle< 2 > &src_tri_local, const Triangle< 3 > &obs_tri, const Triangle< 3 > &src_tri)
Transforms the coordinates of the observation and source triangles into a local coordinate system def...
Definition base.hpp:73
virtual EigMatMN< Complex, obs_num_dof, src_num_dof > compute(const Complex k, const Triangle< 3 > &obs_tri, const Triangle< 3 > &src_tri)=0
Computes operator values for the given observation and source triangles.
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.