16#ifndef dealii_integrators_divergence_h
17#define dealii_integrators_divergence_h
68 const double dx = fe.
JxW(k) * factor;
69 for (
unsigned int i = 0; i < t_dofs; ++i)
72 for (
unsigned int d = 0;
d < dim; ++
d)
73 for (
unsigned int j = 0; j < n_dofs; ++j)
76 M(i, j) +=
dx * du * vv;
91 template <
int dim,
typename number>
96 const double factor = 1.)
106 const double dx = factor * fetest.
JxW(k);
108 for (
unsigned int i = 0; i < t_dofs; ++i)
109 for (
unsigned int d = 0;
d < dim; ++
d)
124 template <
int dim,
typename number>
128 const ArrayView<
const std::vector<double>> &input,
129 const double factor = 1.)
139 const double dx = factor * fetest.
JxW(k);
141 for (
unsigned int i = 0; i < t_dofs; ++i)
142 for (
unsigned int d = 0;
d < dim; ++
d)
172 const double dx = fe.
JxW(k) * factor;
173 for (
unsigned int d = 0;
d < dim; ++
d)
174 for (
unsigned int i = 0; i < t_dofs; ++i)
177 for (
unsigned int j = 0; j < n_dofs; ++j)
180 M(i, j) +=
dx * vv * Du[
d];
195 template <
int dim,
typename number>
200 const double factor = 1.)
210 const double dx = factor * fetest.
JxW(k);
212 for (
unsigned int i = 0; i < t_dofs; ++i)
213 for (
unsigned int d = 0;
d < dim; ++
d)
228 template <
int dim,
typename number>
232 const std::vector<double> &input,
233 const double factor = 1.)
243 const double dx = factor * fetest.
JxW(k);
245 for (
unsigned int i = 0; i < t_dofs; ++i)
246 for (
unsigned int d = 0;
d < dim; ++
d)
275 for (
unsigned int i = 0; i < t_dofs; ++i)
276 for (
unsigned int j = 0; j < n_dofs; ++j)
277 for (
unsigned int d = 0;
d < dim; ++
d)
290 template <
int dim,
typename number>
295 const ArrayView<
const std::vector<double>> &data,
309 for (
unsigned int i = 0; i < t_dofs; ++i)
310 for (
unsigned int d = 0;
d < dim; ++
d)
322 template <
int dim,
typename number>
326 const std::vector<double> &data,
340 for (
unsigned int i = 0; i < t_dofs; ++i)
341 for (
unsigned int d = 0;
d < dim; ++
d)
386 const double dx = factor * fe1.
JxW(k);
387 for (
unsigned int i = 0; i < t_dofs; ++i)
388 for (
unsigned int j = 0; j < n_dofs; ++j)
389 for (
unsigned int d = 0;
d < dim; ++
d)
398 M11(i, j) += .5 *
dx * un1 *
v1;
399 M12(i, j) += .5 *
dx * un2 *
v1;
400 M21(i, j) += .5 *
dx * un1 * v2;
401 M22(i, j) += .5 *
dx * un2 * v2;
440 const double dx = factor * fe1.
JxW(k);
441 for (
unsigned int i = 0; i < n_dofs; ++i)
442 for (
unsigned int j = 0; j < n_dofs; ++j)
443 for (
unsigned int d = 0;
d < dim; ++
d)
454 M11(i, j) +=
dx * un1 * vn1;
455 M12(i, j) +=
dx * un2 * vn1;
456 M21(i, j) +=
dx * un1 * vn2;
457 M22(i, j) +=
dx * un2 * vn2;
481 double div = Du[0][k][0];
482 for (
unsigned int d = 1;
d < dim; ++
d)
484 result += div * div * fe.
JxW(k);
const unsigned int dofs_per_cell
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
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
double JxW(const unsigned int quadrature_point) const
const FiniteElement< dim, spacedim > & get_fe() const
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void u_times_n_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< double > &data, double factor=1.)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
void cell_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< Tensor< 1, dim > > > &input, const double factor=1.)
void u_dot_n_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
void gradient_residual(Vector< number > &result, const FEValuesBase< dim > &fetest, const std::vector< Tensor< 1, dim > > &input, const double factor=1.)
void u_dot_n_residual(Vector< number > &result, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &data, double factor=1.)
void u_dot_n_jump_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double factor=1.)
void gradient_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Library of integrals over cells and faces.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)