15#ifndef dealii_integrators_maxwell_h
16#define dealii_integrators_maxwell_h
101 result[0] = h1[0][1] - h0[1][1];
102 result[1] = h0[0][1] - h1[0][0];
105 result[0] = h1[0][1] + h2[0][2] - h0[1][1] - h0[2][2];
106 result[1] = h2[1][2] + h0[1][0] - h1[2][2] - h1[0][0];
107 result[2] = h0[2][0] + h1[2][1] - h2[0][0] - h2[1][1];
137 result[0] = normal[1] * (g1[0] - g0[1]);
138 result[1] = -normal[0] * (g1[0] - g0[1]);
142 normal[2] * (g2[1] - g0[2]) + normal[1] * (g1[0] - g0[1]);
144 normal[0] * (g0[2] - g1[0]) + normal[2] * (g2[1] - g1[2]);
146 normal[1] * (g1[0] - g2[1]) + normal[0] * (g0[2] - g2[0]);
163 DEAL_II_DEPRECATED_EARLY
void
166 const double factor = 1.)
183 const unsigned int d_max = (dim == 2) ? 1 : dim;
187 const double dx = factor * fe.
JxW(k);
188 for (
unsigned int i = 0; i < n_dofs; ++i)
189 for (
unsigned int j = 0; j < n_dofs; ++j)
190 for (
unsigned int d = 0; d < d_max; ++d)
192 const unsigned int d1 = (d + 1) % dim;
193 const unsigned int d2 = (d + 2) % dim;
200 M(i, j) += dx * cu * cv;
216 DEAL_II_DEPRECATED_EARLY
void
231 const unsigned int d_max = (dim == 2) ? 1 : dim;
235 const double dx = fe.
JxW(k) * factor;
236 for (
unsigned int i = 0; i < t_dofs; ++i)
237 for (
unsigned int j = 0; j < n_dofs; ++j)
238 for (
unsigned int d = 0; d < d_max; ++d)
240 const unsigned int d1 = (d + 1) % dim;
241 const unsigned int d2 = (d + 2) % dim;
246 M(i, j) += dx * cu * vv;
268 const unsigned int face_no,
288 const unsigned int d_max = (dim == 2) ? 1 : dim;
292 const double dx = factor * fe.
JxW(k);
294 for (
unsigned int i = 0; i < n_dofs; ++i)
295 for (
unsigned int j = 0; j < n_dofs; ++j)
299 for (
unsigned int d = 0; d < d_max; ++d)
301 const unsigned int d1 = (d + 1) % dim;
302 const unsigned int d2 = (d + 2) % dim;
315 M(i, j) += dx * (2. * penalty * u * v - cv * u - cu * v);
328 DEAL_II_DEPRECATED_EARLY
void
349 const unsigned int d_max = (dim == 2) ? 1 : dim;
353 const double dx = factor * fe.
JxW(k);
355 for (
unsigned int i = 0; i < n_dofs; ++i)
356 for (
unsigned int j = 0; j < n_dofs; ++j)
357 for (
unsigned int d = 0; d < d_max; ++d)
359 const unsigned int d1 = (d + 1) % dim;
360 const unsigned int d2 = (d + 2) % dim;
367 M(i, j) += dx * u * v;
385 DEAL_II_DEPRECATED_EARLY
inline void
393 const double factor1 = 1.,
394 const double factor2 = -1.)
396 const unsigned int n_dofs = fe1.n_dofs_per_cell();
409 const double nu1 = factor1;
410 const double nu2 = (factor2 < 0) ? factor1 : factor2;
411 const double penalty = .5 * pen * (nu1 + nu2);
423 const unsigned int d_max = (dim == 2) ? 1 : dim;
427 const double dx = fe1.
JxW(k);
429 for (
unsigned int i = 0; i < n_dofs; ++i)
430 for (
unsigned int j = 0; j < n_dofs; ++j)
431 for (
unsigned int d = 0; d < d_max; ++d)
433 const unsigned int d1 = (d + 1) % dim;
434 const unsigned int d2 = (d + 2) % dim;
464 .5 * dx * (2. * penalty * u1 *
v1 - cv1 * u1 - cu1 *
v1);
466 .5 * dx * (2. * penalty *
v1 * u2 - cv1 * u2 - cu2 *
v1);
468 .5 * dx * (2. * penalty * u1 * v2 - cv2 * u1 - cu1 * v2);
470 .5 * dx * (2. * penalty * u2 * v2 - cv2 * u2 - cu2 * v2);
const unsigned int dofs_per_cell
const Tensor< 1, spacedim > & normal_vector(const unsigned int q_point) const
double shape_value_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
const unsigned int n_quadrature_points
Tensor< 1, spacedim > shape_grad_component(const unsigned int i, const unsigned int q_point, const unsigned int component) const
const FiniteElement< dim, spacedim > & get_fe() const
double JxW(const unsigned int q_point) const
unsigned int n_components() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
#define DEAL_II_NOT_IMPLEMENTED()
void nitsche_curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const unsigned int face_no, double penalty, double factor=1.)
void curl_curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
void tangential_trace_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
void curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Tensor< 1, dim > tangential_curl(const Tensor< 1, dim > &g0, const Tensor< 1, dim > &g1, const Tensor< 1, dim > &g2, const Tensor< 1, dim > &normal)
void ip_curl_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const double pen, const double factor1=1., const double factor2=-1.)
Tensor< 1, dim > curl_curl(const Tensor< 2, dim > &h0, const Tensor< 2, dim > &h1, const Tensor< 2, dim > &h2)
Library of integrals over cells and faces.