16 #ifndef dealii__integrators_maxwell_h 17 #define dealii__integrators_maxwell_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/quadrature.h> 23 #include <deal.II/lac/full_matrix.h> 24 #include <deal.II/fe/mapping.h> 25 #include <deal.II/fe/fe_values.h> 26 #include <deal.II/meshworker/dof_info.h> 28 DEAL_II_NAMESPACE_OPEN
102 result[0] = h1[0][1]-h0[1][1];
103 result[1] = h0[0][1]-h1[0][0];
106 result[0] = h1[0][1]+h2[0][2]-h0[1][1]-h0[2][2];
107 result[1] = h2[1][2]+h0[1][0]-h1[2][2]-h1[0][0];
108 result[2] = h0[2][0]+h1[2][1]-h2[0][0]-h2[1][1];
142 result[0] = normal[1] * (g1[0]-g0[1]);
143 result[1] =-normal[0] * (g1[0]-g0[1]);
146 result[0] = normal[2]*(g2[1]-g0[2])+normal[1]*(g1[0]-g0[1]);
147 result[1] = normal[0]*(g0[2]-g1[0])+normal[2]*(g2[1]-g1[2]);
148 result[2] = normal[1]*(g1[0]-g2[1])+normal[0]*(g0[2]-g2[0]);
171 const double factor = 1.)
188 const unsigned int d_max = (dim==2) ? 1 : dim;
192 const double dx = factor * fe.
JxW(k);
193 for (
unsigned int i=0; i<n_dofs; ++i)
194 for (
unsigned int j=0; j<n_dofs; ++j)
195 for (
unsigned int d=0; d<d_max; ++d)
197 const unsigned int d1 = (d+1)%dim;
198 const unsigned int d2 = (d+2)%dim;
203 M(i,j) += dx * cu * cv;
237 const unsigned int d_max = (dim==2) ? 1 : dim;
241 const double dx = fe.
JxW(k) * factor;
242 for (
unsigned int i=0; i<t_dofs; ++i)
243 for (
unsigned int j=0; j<n_dofs; ++j)
244 for (
unsigned int d=0; d<d_max; ++d)
246 const unsigned int d1 = (d+1)%dim;
247 const unsigned int d2 = (d+2)%dim;
251 M(i,j) += dx * cu * vv;
276 const unsigned int face_no,
296 const unsigned int d_max = (dim==2) ? 1 : dim;
300 const double dx = factor * fe.
JxW(k);
302 for (
unsigned int i=0; i<n_dofs; ++i)
303 for (
unsigned int j=0; j<n_dofs; ++j)
304 if (fe.
get_fe().has_support_on_face(i, face_no) && fe.
get_fe().has_support_on_face(j, face_no))
306 for (
unsigned int d=0; d<d_max; ++d)
308 const unsigned int d1 = (d+1)%dim;
309 const unsigned int d2 = (d+2)%dim;
316 M(i,j) += dx*(2.*penalty*u*v - cv*u - cu*v);
353 const unsigned int d_max = (dim==2) ? 1 : dim;
357 const double dx = factor * fe.
JxW(k);
359 for (
unsigned int i=0; i<n_dofs; ++i)
360 for (
unsigned int j=0; j<n_dofs; ++j)
361 for (
unsigned int d=0; d<d_max; ++d)
363 const unsigned int d1 = (d+1)%dim;
364 const unsigned int d2 = (d+2)%dim;
398 const double factor1 = 1.,
399 const double factor2 = -1.)
414 const double nu1 = factor1;
415 const double nu2 = (factor2 < 0) ? factor1 : factor2;
416 const double penalty = .5 * pen * (nu1 + nu2);
428 const unsigned int d_max = (dim==2) ? 1 : dim;
432 const double dx = fe1.
JxW(k);
434 for (
unsigned int i=0; i<n_dofs; ++i)
435 for (
unsigned int j=0; j<n_dofs; ++j)
436 for (
unsigned int d=0; d<d_max; ++d)
438 const unsigned int d1 = (d+1)%dim;
439 const unsigned int d2 = (d+2)%dim;
452 M11(i,j) += .5*dx*(2.*penalty*u1*v1 - cv1*u1 - cu1*v1);
453 M12(i,j) += .5*dx*(2.*penalty*v1*u2 - cv1*u2 - cu2*v1);
454 M21(i,j) += .5*dx*(2.*penalty*u1*v2 - cv2*u1 - cu1*v2);
455 M22(i,j) += .5*dx*(2.*penalty*u2*v2 - cv2*u2 - cu2*v2);
465 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
const unsigned int dofs_per_cell
const FiniteElement< dim, spacedim > & get_fe() const
Tensor< 1, dim > curl_curl(const Tensor< 2, dim > &h0, const Tensor< 2, dim > &h1, const Tensor< 2, dim > &h2)
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 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
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()
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const