16 #ifndef dealii_integrators_maxwell_h 17 #define dealii_integrators_maxwell_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/quadrature.h> 25 #include <deal.II/fe/fe_values.h> 26 #include <deal.II/fe/mapping.h> 28 #include <deal.II/lac/full_matrix.h> 30 #include <deal.II/meshworker/dof_info.h> 32 DEAL_II_NAMESPACE_OPEN
105 result[0] = h1[0][1] - h0[1][1];
106 result[1] = h0[0][1] - h1[0][0];
109 result[0] = h1[0][1] + h2[0][2] - h0[1][1] - h0[2][2];
110 result[1] = h2[1][2] + h0[1][0] - h1[2][2] - h1[0][0];
111 result[2] = h0[2][0] + h1[2][1] - h2[0][0] - h2[1][1];
144 result[0] = normal[1] * (g1[0] - g0[1]);
145 result[1] = -normal[0] * (g1[0] - g0[1]);
149 normal[2] * (g2[1] - g0[2]) + normal[1] * (g1[0] - g0[1]);
151 normal[0] * (g0[2] - g1[0]) + normal[2] * (g2[1] - g1[2]);
153 normal[1] * (g1[0] - g2[1]) + normal[0] * (g0[2] - g2[0]);
176 const double factor = 1.)
193 const unsigned int d_max = (dim == 2) ? 1 : dim;
197 const double dx = factor * fe.
JxW(k);
198 for (
unsigned int i = 0; i < n_dofs; ++i)
199 for (
unsigned int j = 0; j < n_dofs; ++j)
200 for (
unsigned int d = 0; d < d_max; ++d)
202 const unsigned int d1 = (d + 1) % dim;
203 const unsigned int d2 = (d + 2) % dim;
210 M(i, j) += dx * cu * cv;
244 const unsigned int d_max = (dim == 2) ? 1 : dim;
248 const double dx = fe.
JxW(k) * factor;
249 for (
unsigned int i = 0; i < t_dofs; ++i)
250 for (
unsigned int j = 0; j < n_dofs; ++j)
251 for (
unsigned int d = 0; d < d_max; ++d)
253 const unsigned int d1 = (d + 1) % dim;
254 const unsigned int d2 = (d + 2) % dim;
259 M(i, j) += dx * cu * vv;
284 const unsigned int face_no,
304 const unsigned int d_max = (dim == 2) ? 1 : dim;
308 const double dx = factor * fe.
JxW(k);
310 for (
unsigned int i = 0; i < n_dofs; ++i)
311 for (
unsigned int j = 0; j < n_dofs; ++j)
312 if (fe.
get_fe().has_support_on_face(i, face_no) &&
313 fe.
get_fe().has_support_on_face(j, face_no))
315 for (
unsigned int d = 0; d < d_max; ++d)
317 const unsigned int d1 = (d + 1) % dim;
318 const unsigned int d2 = (d + 2) % dim;
331 M(i, j) += dx * (2. * penalty * u * v - cv * u - cu * v);
368 const unsigned int d_max = (dim == 2) ? 1 : dim;
372 const double dx = factor * fe.
JxW(k);
374 for (
unsigned int i = 0; i < n_dofs; ++i)
375 for (
unsigned int j = 0; j < n_dofs; ++j)
376 for (
unsigned int d = 0; d < d_max; ++d)
378 const unsigned int d1 = (d + 1) % dim;
379 const unsigned int d2 = (d + 2) % dim;
386 M(i, j) += dx * u * v;
415 const double factor1 = 1.,
416 const double factor2 = -1.)
431 const double nu1 = factor1;
432 const double nu2 = (factor2 < 0) ? factor1 : factor2;
433 const double penalty = .5 * pen * (nu1 + nu2);
445 const unsigned int d_max = (dim == 2) ? 1 : dim;
449 const double dx = fe1.
JxW(k);
451 for (
unsigned int i = 0; i < n_dofs; ++i)
452 for (
unsigned int j = 0; j < n_dofs; ++j)
453 for (
unsigned int d = 0; d < d_max; ++d)
455 const unsigned int d1 = (d + 1) % dim;
456 const unsigned int d2 = (d + 2) % dim;
486 .5 * dx * (2. * penalty * u1 * v1 - cv1 * u1 - cu1 * v1);
488 .5 * dx * (2. * penalty * v1 * u2 - cv1 * u2 - cu2 * v1);
490 .5 * dx * (2. * penalty * u1 * v2 - cv2 * u1 - cu1 * v2);
492 .5 * dx * (2. * penalty * u2 * v2 - cv2 * u2 - cu2 * v2);
502 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
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)
const unsigned int dofs_per_cell
const FiniteElement< dim, spacedim > & get_fe() const
void tangential_trace_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
Library of integrals over cells and faces.
void curl_curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
#define Assert(cond, exc)
void curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
const unsigned int n_quadrature_points
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
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.)
void nitsche_curl_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const unsigned int face_no, double penalty, double factor=1.)
double JxW(const unsigned int quadrature_point) const
static ::ExceptionBase & ExcNotImplemented()
Tensor< 1, dim > curl_curl(const Tensor< 2, dim > &h0, const Tensor< 2, dim > &h1, const Tensor< 2, dim > &h2)
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const