16 #ifndef dealii_integrators_l2_h 17 #define dealii_integrators_l2_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
56 const double factor = 1.)
59 const unsigned int n_components = fe.
get_fe().n_components();
63 const double dx = fe.
JxW(k) * factor;
64 for (
unsigned int i=0; i<n_dofs; ++i)
67 for (
unsigned int d=0; d<n_components; ++d)
74 for (
unsigned int j=i+1; j<n_dofs; ++j)
77 for (
unsigned int d=0; d<n_components; ++d)
108 const std::vector<double> &weights)
111 const unsigned int n_components = fe.
get_fe().n_components();
118 const double dx = fe.
JxW(k) * weights[k];
119 for (
unsigned int i=0; i<n_dofs; ++i)
122 for (
unsigned int d=0; d<n_components; ++d)
129 for (
unsigned int j=i+1; j<n_dofs; ++j)
132 for (
unsigned int d=0; d<n_components; ++d)
152 template <
int dim,
typename number>
156 const std::vector<double> &input,
157 const double factor = 1.)
165 for (
unsigned int i=0; i<n_dofs; ++i)
177 template <
int dim,
typename number>
181 const VectorSlice<
const std::vector<std::vector<double> > > &input,
182 const double factor = 1.)
185 const unsigned int n_components = input.size();
191 for (
unsigned int i=0; i<n_dofs; ++i)
192 for (
unsigned int d=0; d<n_components; ++d)
216 const double factor1 = 1.,
217 const double factor2 = 1.)
221 const unsigned int n_components = fe1.
get_fe().n_components();
237 const double dx = fe1.
JxW(k);
239 for (
unsigned int i=0; i<n1_dofs; ++i)
240 for (
unsigned int j=0; j<n1_dofs; ++j)
241 for (
unsigned int d=0; d<n_components; ++d)
248 M11(i,j) += dx * u1*v1;
249 M12(i,j) += dx * u2*v1;
250 M21(i,j) += dx * u1*v2;
251 M22(i,j) += dx * u2*v2;
258 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
const unsigned int dofs_per_cell
const FiniteElement< dim, spacedim > & get_fe() const
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.)
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 weighted_mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const std::vector< double > &weights)
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
#define Assert(cond, exc)
const unsigned int n_quadrature_points
double JxW(const unsigned int quadrature_point) const
static ::ExceptionBase & ExcNotImplemented()
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)