16#ifndef dealii_integrators_grad_div_h
17#define dealii_integrators_grad_div_h
64 const double dx = factor * fe.
JxW(k);
65 for (
unsigned int i = 0; i < n_dofs; ++i)
66 for (
unsigned int j = 0; j < n_dofs; ++j)
73 M(i, j) += dx * divu * divv;
84 template <
int dim,
typename number>
89 const double factor = 1.)
98 const double dx = factor * fetest.
JxW(k);
99 for (
unsigned int i = 0; i < n_dofs; ++i)
104 for (
unsigned int d = 0; d < dim; ++d)
105 du += input[d][k][d];
107 result(i) += dx * du * divv;
135 const double dx = factor * fe.
JxW(k);
137 for (
unsigned int i = 0; i < n_dofs; ++i)
138 for (
unsigned int j = 0; j < n_dofs; ++j)
144 double un = 0., vn = 0.;
145 for (
unsigned int d = 0; d < dim; ++d)
151 M(i, j) += dx * 2. * penalty * un * vn;
152 M(i, j) -= dx * (divu * vn + divv * un);
176 const ArrayView<
const std::vector<double>> & input,
178 const ArrayView<
const std::vector<double>> & data,
190 const double dx = factor * fe.
JxW(k);
195 for (
unsigned int d = 0; d < dim; ++d)
197 umgn += (input[d][k] - data[d][k]) * n[d];
198 divu += Dinput[d][k][d];
201 for (
unsigned int i = 0; i < n_dofs; ++i)
206 for (
unsigned int d = 0; d < dim; ++d)
210 dx * (2. * penalty * umgn * vn - divv * umgn - divu * vn);
230 double factor2 = -1.)
242 const double fi = factor1;
243 const double fe = (factor2 < 0) ? factor1 : factor2;
244 const double f = .5 * (fi + fe);
248 const double dx = fe1.
JxW(k);
250 for (
unsigned int i = 0; i < n_dofs; ++i)
251 for (
unsigned int j = 0; j < n_dofs; ++j)
266 for (
unsigned int d = 0; d < dim; ++d)
274 dx * (-.5 * fi * divvi * uni - .5 * fi * divui * vni +
275 f * penalty * uni * vni);
277 dx * (.5 * fi * divvi * une - .5 * fe * divue * vni -
278 f * penalty * vni * une);
280 dx * (-.5 * fe * divve * uni + .5 * fi * divui * vne -
281 f * penalty * uni * vne);
283 dx * (.5 * fe * divve * une + .5 * fe * divue * vne +
284 f * penalty * une * vne);
306 const ArrayView<
const std::vector<double>> & input1,
308 const ArrayView<
const std::vector<double>> & input2,
311 double int_factor = 1.,
312 double ext_factor = -1.)
322 const double fi = int_factor;
323 const double fe = (ext_factor < 0) ? int_factor : ext_factor;
324 const double penalty = .5 * pen * (fi + fe);
329 const double dx = fe1.
JxW(k);
335 for (
unsigned int d = 0; d < dim; ++d)
337 uni += input1[d][k] * n[d];
338 une += input2[d][k] * n[d];
339 divui += Dinput1[d][k][d];
340 divue += Dinput2[d][k][d];
343 for (
unsigned int i = 0; i < n1; ++i)
351 for (
unsigned int d = 0; d < dim; ++d)
357 result1(i) += dx * (-.5 * fi * divvi * uni -
358 .5 * fi * divui * vni + penalty * uni * vni);
359 result1(i) += dx * (.5 * fi * divvi * une -
360 .5 * fe * divue * vni - penalty * vni * une);
361 result2(i) += dx * (-.5 * fe * divve * uni +
362 .5 * fi * divui * vne - penalty * uni * vne);
363 result2(i) += dx * (.5 * fe * divve * une +
364 .5 * fe * divue * vne + penalty * une * vne);
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
const FiniteElement< dim, spacedim > & get_fe() const
double JxW(const unsigned int q_point) const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
#define AssertDimension(dim1, dim2)
void nitsche_residual(Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double > > &input, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput, const ArrayView< const std::vector< double > > &data, double penalty, double factor=1.)
void nitsche_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, double factor=1.)
void ip_residual(Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const ArrayView< const std::vector< double > > &input1, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput1, const ArrayView< const std::vector< double > > &input2, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
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 cell_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< Tensor< 1, dim > > > &input, const double factor=1.)
Library of integrals over cells and faces.