|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_integrators_l2_h
17 #define dealii_integrators_l2_h
65 const double factor = 1.)
72 const double dx = fe.
JxW(k) * factor;
73 for (
unsigned int i = 0; i < n_dofs; ++i)
82 for (
unsigned int j = i + 1; j < n_dofs; ++j)
119 const std::vector<double> &weights)
129 const double dx = fe.
JxW(k) * weights[k];
130 for (
unsigned int i = 0; i < n_dofs; ++i)
139 for (
unsigned int j = i + 1; j < n_dofs; ++j)
169 template <
int dim,
typename number>
173 const std::vector<double> &input,
174 const double factor = 1.)
182 for (
unsigned int i = 0; i < n_dofs; ++i)
183 result(i) += fe.
JxW(k) * factor * input[k] * fe.
shape_value(i, k);
202 template <
int dim,
typename number>
206 const ArrayView<
const std::vector<double>> &input,
207 const double factor = 1.)
216 for (
unsigned int i = 0; i < n_dofs; ++i)
218 result(i) += fe.
JxW(k) * factor *
261 const double factor1 = 1.,
262 const double factor2 = 1.)
282 const double dx = fe1.
JxW(k);
284 for (
unsigned int i = 0; i < n1_dofs; ++i)
285 for (
unsigned int j = 0; j < n1_dofs; ++j)
297 M11(i, j) +=
dx * u1 * v1;
298 M12(i, j) +=
dx * u2 * v1;
299 M21(i, j) +=
dx * u1 * v2;
300 M22(i, j) +=
dx * u2 * v2;
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.)
static ::ExceptionBase & ExcNotImplemented()
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double >> &input, const double factor=1.)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const unsigned int dofs_per_cell
#define DEAL_II_NAMESPACE_OPEN
#define AssertDimension(dim1, dim2)
double JxW(const unsigned int quadrature_point) const
Library of integrals over cells and faces.
void L2(Vector< number > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const double factor=1.)
#define Assert(cond, exc)
const FiniteElement< dim, spacedim > & get_fe() const
const unsigned int n_quadrature_points
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)
#define DEAL_II_NAMESPACE_CLOSE