16#ifndef dealii_integrators_l2_h
17#define dealii_integrators_l2_h
60 const double factor = 1.)
67 const double dx = fe.
JxW(k) * factor;
68 for (
unsigned int i = 0; i < n_dofs; ++i)
71 for (
unsigned int d = 0; d < n_components; ++d)
77 for (
unsigned int j = i + 1; j < n_dofs; ++j)
80 for (
unsigned int d = 0; d < n_components; ++d)
111 const std::vector<double> &weights)
121 const double dx = fe.
JxW(k) * weights[k];
122 for (
unsigned int i = 0; i < n_dofs; ++i)
125 for (
unsigned int d = 0; d < n_components; ++d)
131 for (
unsigned int j = i + 1; j < n_dofs; ++j)
134 for (
unsigned int d = 0; d < n_components; ++d)
158 template <
int dim,
typename number>
162 const std::vector<double> &input,
163 const double factor = 1.)
171 for (
unsigned int i = 0; i < n_dofs; ++i)
172 result(i) += fe.
JxW(k) * factor * input[k] * fe.
shape_value(i, k);
188 template <
int dim,
typename number>
192 const ArrayView<
const std::vector<double>> &input,
193 const double factor = 1.)
196 const unsigned int n_components = input.size();
202 for (
unsigned int i = 0; i < n_dofs; ++i)
203 for (
unsigned int d = 0; d < n_components; ++d)
204 result(i) += fe.
JxW(k) * factor *
244 const double factor1 = 1.,
245 const double factor2 = 1.)
247 const unsigned int n1_dofs = fe1.n_dofs_per_cell();
248 const unsigned int n2_dofs = fe2.n_dofs_per_cell();
265 const double dx = fe1.
JxW(k);
267 for (
unsigned int i = 0; i < n1_dofs; ++i)
268 for (
unsigned int j = 0; j < n1_dofs; ++j)
269 for (
unsigned int d = 0; d < n_components; ++d)
280 M11(i, j) += dx * u1 *
v1;
281 M12(i, j) += dx * u2 *
v1;
282 M21(i, j) += dx * u1 * v2;
283 M22(i, j) += dx * u2 * v2;
const unsigned int dofs_per_cell
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
const double & shape_value(const unsigned int i, const unsigned int q_point) const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
void jump_matrix(FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const double factor1=1., const double factor2=1.)
void weighted_mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &weights)
Library of integrals over cells and faces.