From: Marco Feder Date: Tue, 23 Nov 2021 19:43:18 +0000 (+0100) Subject: step-22 now uses extractors in local rhs X-Git-Tag: v9.4.0-rc1~804^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F12988%2Fhead;p=dealii.git step-22 now uses extractors in local rhs --- diff --git a/doc/news/changes/minor/20211123MarcoFeder b/doc/news/changes/minor/20211123MarcoFeder new file mode 100644 index 0000000000..275b07482c --- /dev/null +++ b/doc/news/changes/minor/20211123MarcoFeder @@ -0,0 +1,3 @@ +Fixed: The step-22 tutorial now uses FEValuesExtractors also for the right-hand side +
+(Marco Feder, 2021/11/23) \ No newline at end of file diff --git a/examples/step-22/step-22.cc b/examples/step-22/step-22.cc index 16ccdf0015..1ab4c866c7 100644 --- a/examples/step-22/step-22.cc +++ b/examples/step-22/step-22.cc @@ -25,6 +25,7 @@ #include #include #include +#include #include #include @@ -42,6 +43,7 @@ #include #include + #include #include #include @@ -216,37 +218,43 @@ namespace Step22 // We implement similar functions for the right hand side which for the - // current example is simply zero: + // current example is simply zero. In the current case, what we want to + // represent is the right hand side of the velocity equation, that is, the + // body forces acting on the fluid. (The right hand side of the divergence + // equation is always zero, since we are dealing with an incompressible + // medium.) As such, the right hand side is represented by a `Tensor<1,dim>` + // object: a dim-dimensional vector, and the right tool to implement such a + // function is to derive it from `TensorFunction<1,dim>`: template - class RightHandSide : public Function + class RightHandSide : public TensorFunction<1, dim> { public: RightHandSide() - : Function(dim + 1) + : TensorFunction<1, dim>() {} - virtual double value(const Point & p, - const unsigned int component = 0) const override; + virtual Tensor<1, dim> value(const Point &p) const override; - virtual void vector_value(const Point &p, - Vector & value) const override; + virtual void value_list(const std::vector> &p, + std::vector> &value) const override; }; template - double RightHandSide::value(const Point & /*p*/, - const unsigned int /*component*/) const + Tensor<1, dim> RightHandSide::value(const Point & /*p*/) const { - return 0; + return Tensor<1, dim>(); } template - void RightHandSide::vector_value(const Point &p, - Vector & values) const + void RightHandSide::value_list(const std::vector> &vp, + std::vector> &values) const { - for (unsigned int c = 0; c < this->n_components; ++c) - values(c) = RightHandSide::value(p, c); + for (unsigned int c = 0; c < vp.size(); ++c) + { + values[c] = RightHandSide::value(vp[c]); + } } @@ -629,7 +637,7 @@ namespace Step22 std::vector local_dof_indices(dofs_per_cell); const RightHandSide right_hand_side; - std::vector> rhs_values(n_q_points, Vector(dim + 1)); + std::vector> rhs_values(n_q_points, Tensor<1, dim>()); // Next, we need two objects that work as extractors for the FEValues // object. Their use is explained in detail in the report on @ref @@ -661,6 +669,7 @@ namespace Step22 // the outer loop. std::vector> symgrad_phi_u(dofs_per_cell); std::vector div_phi_u(dofs_per_cell); + std::vector> phi_u(dofs_per_cell); std::vector phi_p(dofs_per_cell); for (const auto &cell : dof_handler.active_cell_iterators()) @@ -670,8 +679,8 @@ namespace Step22 local_preconditioner_matrix = 0; local_rhs = 0; - right_hand_side.vector_value_list(fe_values.get_quadrature_points(), - rhs_values); + right_hand_side.value_list(fe_values.get_quadrature_points(), + rhs_values); for (unsigned int q = 0; q < n_q_points; ++q) { @@ -680,6 +689,7 @@ namespace Step22 symgrad_phi_u[k] = fe_values[velocities].symmetric_gradient(k, q); div_phi_u[k] = fe_values[velocities].divergence(k, q); + phi_u[k] = fe_values[velocities].value(k, q); phi_p[k] = fe_values[pressure].value(k, q); } @@ -729,22 +739,17 @@ namespace Step22 // is overloaded for symmetric tensors, yielding the scalar // product between the two tensors. // - // For the right-hand side we use the fact that the shape - // functions are only non-zero in one component (because our - // elements are primitive). Instead of multiplying the tensor - // representing the dim+1 values of shape function i with the - // whole right-hand side vector, we only look at the only - // non-zero component. The function - // FiniteElement::system_to_component_index will return - // which component this shape function lives in (0=x velocity, - // 1=y velocity, 2=pressure in 2d), which we use to pick out - // the correct component of the right-hand side vector to - // multiply with. - const unsigned int component_i = - fe.system_to_component_index(i).first; - local_rhs(i) += (fe_values.shape_value(i, q) // (phi_u_i(x_q) - * rhs_values[q](component_i)) // * f(x_q)) - * fe_values.JxW(q); // * dx + // For the right-hand side, we need to multiply the (vector of) + // velocity shape functions with the vector of body force + // right-hand side components, both evaluated at the current + // quadrature point. We have implemented the body forces as a + // `TensorFunction<1,dim>`, so its values at quadrature points + // are already tensors for which the application of `operator*` + // against the velocity components of the shape function results + // in the dot product, as intended. + local_rhs(i) += phi_u[i] // phi_u_i(x_q) + * rhs_values[q] // * f(x_q) + * fe_values.JxW(q); // * dx } }