16 #ifndef dealii__integrators_laplace_h 17 #define dealii__integrators_laplace_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
55 const double factor = 1.)
58 const unsigned int n_components = fe.
get_fe().n_components();
62 const double dx = fe.
JxW(k) * factor;
63 for (
unsigned int i=0; i<n_dofs; ++i)
66 for (
unsigned int d=0; d<n_components; ++d)
72 for (
unsigned int j=i+1; j<n_dofs; ++j)
75 for (
unsigned int d=0; d<n_components; ++d)
105 for (
unsigned int k=0; k<nq; ++k)
107 const double dx = factor * fe.
JxW(k);
108 for (
unsigned int i=0; i<n_dofs; ++i)
109 result(i) += dx * (input[k] * fe.
shape_grad(i,k));
129 const unsigned int n_comp = fe.
get_fe().n_components();
134 for (
unsigned int k=0; k<nq; ++k)
136 const double dx = factor * fe.
JxW(k);
137 for (
unsigned int i=0; i<n_dofs; ++i)
138 for (
unsigned int d=0; d<n_comp; ++d)
168 const unsigned int n_comp = fe.
get_fe().n_components();
175 const double dx = fe.
JxW(k) * factor;
177 for (
unsigned int i=0; i<n_dofs; ++i)
178 for (
unsigned int j=0; j<n_dofs; ++j)
179 for (
unsigned int d=0; d<n_comp; ++d)
215 const double dx = fe.
JxW(k) * factor;
217 for (
unsigned int i=0; i<n_dofs; ++i)
218 for (
unsigned int j=0; j<n_dofs; ++j)
225 for (
unsigned int d=0; d<dim; ++d)
233 for (
unsigned int d=0; d<dim; ++d)
240 M(i,j) += dx *(2.*penalty*u_t*v_t-dnu_t*v_t-dnv_t*u_t);
265 const std::vector<double> &input,
267 const std::vector<double> &data,
278 const double dx = factor * fe.
JxW(k);
280 for (
unsigned int i=0; i<n_dofs; ++i)
283 const double dnu = Dinput[k] * n;
285 const double u= input[k];
286 const double g= data[k];
288 result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
314 const VectorSlice<
const std::vector<std::vector<double> > > &input,
316 const VectorSlice<
const std::vector<std::vector<double> > > &data,
321 const unsigned int n_comp = fe.
get_fe().n_components();
328 const double dx = factor * fe.
JxW(k);
330 for (
unsigned int i=0; i<n_dofs; ++i)
331 for (
unsigned int d=0; d<n_comp; ++d)
334 const double dnu = Dinput[d][k] * n;
336 const double u= input[d][k];
337 const double g= data[d][k];
339 result(i) += dx*(2.*penalty*(u-g)*v - dnv*(u-g) - dnu*v);
373 double factor2 = -1.)
385 const double nui = factor1;
386 const double nue = (factor2 < 0) ? factor1 : factor2;
387 const double nu = .5*(nui+nue);
391 const double dx = fe1.
JxW(k);
393 for (
unsigned int d=0; d<fe1.
get_fe().n_components(); ++d)
395 for (
unsigned int i=0; i<n_dofs; ++i)
397 for (
unsigned int j=0; j<n_dofs; ++j)
407 M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
408 M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
409 M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
410 M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
441 double factor2 = -1.)
455 const double nui = factor1;
456 const double nue = (factor2 < 0) ? factor1 : factor2;
457 const double nu = .5*(nui+nue);
461 const double dx = fe1.
JxW(k);
463 for (
unsigned int i=0; i<n_dofs; ++i)
468 for (
unsigned int j=0; j<n_dofs; ++j)
475 double ngradu1n = 0.;
476 double ngradv1n = 0.;
477 double ngradu2n = 0.;
478 double ngradv2n = 0.;
480 for (
unsigned int d=0; d<dim; ++d)
496 for (
unsigned int d=0; d<dim; ++d)
510 M11(i,j) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+nu*penalty*ui*vi);
511 M12(i,j) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-nu*penalty*vi*ue);
512 M21(i,j) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-nu*penalty*ui*ve);
513 M22(i,j) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+nu*penalty*ue*ve);
537 const std::vector<double> &input1,
539 const std::vector<double> &input2,
542 double int_factor = 1.,
543 double ext_factor = -1.)
550 const double nui = int_factor;
551 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
552 const double penalty = .5 * pen * (nui + nue);
558 const double dx = fe1.
JxW(k);
561 for (
unsigned int i=0; i<n_dofs; ++i)
565 const double dnvi = Dvi * n;
568 const double dnve = Dve * n;
570 const double ui = input1[k];
572 const double dnui = Dui * n;
573 const double ue = input2[k];
575 const double dnue = Due * n;
577 result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
578 result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
579 result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
580 result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
604 const VectorSlice<
const std::vector<std::vector<double> > > &input1,
606 const VectorSlice<
const std::vector<std::vector<double> > > &input2,
609 double int_factor = 1.,
610 double ext_factor = -1.)
612 const unsigned int n_comp = fe1.
get_fe().n_components();
620 const double nui = int_factor;
621 const double nue = (ext_factor < 0) ? int_factor : ext_factor;
622 const double penalty = .5 * pen * (nui + nue);
627 const double dx = fe1.
JxW(k);
630 for (
unsigned int i=0; i<n1; ++i)
631 for (
unsigned int d=0; d<n_comp; ++d)
635 const double dnvi = Dvi * n;
638 const double dnve = Dve * n;
640 const double ui = input1[d][k];
642 const double dnui = Dui * n;
643 const double ue = input2[d][k];
645 const double dnue = Due * n;
647 result1(i) += dx*(-.5*nui*dnvi*ui-.5*nui*dnui*vi+penalty*ui*vi);
648 result1(i) += dx*( .5*nui*dnvi*ue-.5*nue*dnue*vi-penalty*vi*ue);
649 result2(i) += dx*(-.5*nue*dnve*ui+.5*nui*dnui*ve-penalty*ui*ve);
650 result2(i) += dx*( .5*nue*dnve*ue+.5*nue*dnue*ve+penalty*ue*ve);
671 template <
int dim,
int spacedim,
typename number>
680 const unsigned int deg1sq = (deg1 == 0) ? 1 : deg1 * (deg1+1);
681 const unsigned int deg2sq = (deg2 == 0) ? 1 : deg2 * (deg2+1);
683 double penalty1 = deg1sq / dinfo1.
cell->extent_in_direction(normal1);
684 double penalty2 = deg2sq / dinfo2.
cell->extent_in_direction(normal2);
685 if (dinfo1.
cell->has_children() ^ dinfo2.
cell->has_children())
691 const double penalty = 0.5*(penalty1 + penalty2);
698 DEAL_II_NAMESPACE_CLOSE
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
#define AssertDimension(dim1, dim2)
void ip_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
void nitsche_tangential_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
const unsigned int dofs_per_cell
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
#define AssertVectorVectorDimension(vec, dim1, dim2)
const FiniteElement< dim, spacedim > & get_fe() const
void cell_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, 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 nitsche_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const std::vector< Tensor< 1, dim > > &Dinput, const std::vector< double > &data, double penalty, double factor=1.)
void ip_residual(Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const std::vector< double > &input1, const std::vector< Tensor< 1, dim > > &Dinput1, const std::vector< double > &input2, const std::vector< Tensor< 1, dim > > &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const unsigned int n_quadrature_points
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
double compute_penalty(const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo1, const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo2, unsigned int deg1, unsigned int deg2)
double JxW(const unsigned int quadrature_point) const
Triangulation< dim, spacedim >::face_iterator face
The current face.
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
void ip_tangential_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
static ::ExceptionBase & ExcInternalError()
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