From 564581da1e247802514430b4c92dc776abe02f58 Mon Sep 17 00:00:00 2001 From: bangerth Date: Sun, 14 Dec 2008 04:50:01 +0000 Subject: [PATCH] Document new FEValuesViews features. git-svn-id: https://svn.dealii.org/trunk@17944 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/doc/doxygen/headers/vector_valued.h | 151 ++++++++++++++++++-- deal.II/doc/news/changes.h | 16 ++- 2 files changed, 156 insertions(+), 11 deletions(-) diff --git a/deal.II/doc/doxygen/headers/vector_valued.h b/deal.II/doc/doxygen/headers/vector_valued.h index 4fbda7dc07..08a9800db6 100644 --- a/deal.II/doc/doxygen/headers/vector_valued.h +++ b/deal.II/doc/doxygen/headers/vector_valued.h @@ -47,6 +47,7 @@ *
  • @ref VVAssembling "Assembling linear systems" *
  • @ref VVAlternative "An alternative approach" *
  • @ref VVBlockSolvers "Block solvers" + *
  • @ref VVExtracting "Extracting data from solutions" * * * @@ -223,16 +224,23 @@ * copy of the element used for the pressure $p$. * * In this manner, we can combine the whole vector-valued element from its - * individual components. On a sidenote it may also be worth mentioning why we - * want to also mirror our logical construction of the weak form in our - * program: by creating a finite element object as shown above and then - * handing it off to a DoFHandler object, we get an enumeration of all - * degrees of freedom at once, whether velocity or pressure. Initially, they - * are enumerated in a rather random order, but they can be arranged such that - * first all velocity components are enumerated and then all pressures by, for - * example, calling DoFRenumbering::component_wise . By and large, many things - * become much simpler if we can rely on the fact that there is only a single - * DoFHandler that describes everything. + * individual components. However, it is worth pointing out that whichever + * way you construct your finite element object, this object doesn't really + * know what it represents. For example, for + * @code + * FESystem finite_element (FE_Q(1), dim, + * FE_Q(1), 1); + * @endcode + * the constructed object really only knows that it has dim+1 + * vector components. It has no notion, however, which of these components + * represent scalar fields (e.g. temperature, pressure, concentration) and/or + * if any of its components are parts of vector fields (velocities, + * displacements) or tensors (e.g. stresses). As a consequence, the FEValues + * objects we use below to evaluate finite element shape functions at + * quadrature points only knows that it has a finite element with a number + * of vector components, but doesn't know how to group them. We will show + * how to give these components a logical connection using the + * FEValuesExtractors classes. * * * @anchor VVAssembling @@ -738,6 +746,129 @@ scalar_product (const Tensor<2,dim> &u, * therefore allows us to quite easily write rather complicated solvers making * use of the block structure of a linear system. * + * + * + * @anchor VVExtracting + *

    Extracting data from solutions

    + * + * Once one has computed a solution, it is often necessary to evaluate it at + * quadrature points, for example to evaluate nonlinear residuals for the next + * Newton iteration, to evaluate the finite element residual for error + * estimators, or to compute the right hand side for the next time step in + * a time dependent problem. + * + * The way this is done us to again use an FEValues object to evaluate + * the shape functions at quadrature points, and with those also the + * values of a finite element function. For the example of the mixed + * Laplace problem above, consider the following code after solving: + * @code + std::vector > local_solution_values (n_q_points, + Vector (dim+1)); + + typename DoFHandler::active_cell_iterator + cell = dof_handler.begin_active(), + endc = dof_handler.end(); + for (; cell!=endc; ++cell) + { + fe_values.reinit (cell); + + fe_values.get_function_values (solution, + local_solution_values); + * @endcode + * + * After this, the variable local_solution_values is a list + * of vectors of a length equal to the number of quadrature points we + * have initialized the FEValues object with; each of these vectors + * has dim+1 elements containing the values of the + * dim velocities and the one pressure at a quadrature point. + * + * We can use these values to then construct other things like residuals. + * However, the construct is a bit awkward. First, we have a + * std::vector of dealii::Vectors, which always + * looks strange. It is also inefficient because it implies dynamic memory + * allocation for the outer vector as well as for all the inner vectors. + * Secondly, maybe we are only interested in the velocities, + * for example to solve an advection problem in a second stage (as, for + * example, in @ref step_21 "step-21" or @ref step_31 "step-31"). In that + * case, one would have to hand-extract these values like so: + * @code + for (unsigned int q=0; q velocity; + for (unsigned int d=0; dTensor@<1,dim@> because the + * velocity is a quantity characterized by dim elements that + * have certain transformation properties under rotations of the coordinate + * system. + * + * This code can be written more elegantly and efficiently using code like + * the following: + * @code + std::vector > local_velocity_values (n_q_points); + + const FEValuesExtractors::Vector velocities (0); + + typename DoFHandler::active_cell_iterator + cell = dof_handler.begin_active(), + endc = dof_handler.end(); + for (; cell!=endc; ++cell) + { + fe_values.reinit (cell); + + fe_values[velocities].get_function_values (solution, + local_velocity_values); + * @endcode + * + * As a result, we here get the velocities right away, and in the + * right data type (because we have described, using the extractor, + * that the first dim components of the finite element + * belong together, forming a tensor). The code is also more efficient: + * it requires less dynamic memory allocation because the Tensor + * class allocates its components as member variables rather than on + * the heap, and we save cycles because we don't even bother computing + * the values of the pressure variable at the quadrature points. On + * the other hand, if we had been interested in only the pressure and + * not the velocities, then the following code extracting scalar + * values would have done: + * @code + std::vector local_pressure_values (n_q_points); + + const FEValuesExtractors::Scalar pressure (dim); + + typename DoFHandler::active_cell_iterator + cell = dof_handler.begin_active(), + endc = dof_handler.end(); + for (; cell!=endc; ++cell) + { + fe_values.reinit (cell); + + fe_values[pressure].get_function_values (solution, + local_pressure_values); + * @endcode + * + * In similar cases, one sometimes needs the gradients or second + * derivatives of the solution, or of individual scalar or vector + * components. To get at those of all components of the solution, + * the functions FEValuesBase::get_function_gradients and + * FEValuesBase::get_function_hessians are the equivalent of the + * function FEValuesBase::get_function_values used above. + * + * Likewise, to extract the gradients of scalar components, + * FEValuesViews::Scalar::get_function_gradients and + * FEValuesViews::Scalar::get_function_hessians do the job. + * For vector-(tensor-)valued quantities, there are functions + * FEValuesViews::Vector::get_function_gradients and + * FEValuesViews::Vector::get_function_hessians, and in + * addition + * FEValuesViews::Vector::get_function_symmetric_gradients and + * FEValuesViews::Vector::get_function_divergences. + * + * * @ingroup feall feaccess */ diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index ed988f0d1b..70a46561be 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -598,6 +598,20 @@ inconvenience this causes.

    deal.II

      +
    1. +

      + Improved: The FEValuesViews objects that one gets when writing things + like fe_values[velocities] have acquired the ability to + extract individual components of solutions at quadrature points by + writing code like + fe_values[velocities].get_function_values (global_solution, + local_velocity_values). See the new section at the end of the + @ref vector_valued module for an overview of this facility. The + @ref step_31 "step-31" program has also been updated to use this. +
      + (WB 2008/12/13) +

      +
    2. Improved: The FEValuesViews objects that one gets when writing things @@ -611,7 +625,7 @@ inconvenience this causes. pre-generated object. This turns an $O(N)$ effort into an $O(1)$ effort, where $N$ is the number of cells.
      - (Luca Heltai 2008/12/08) + (WB 2008/12/10)

    3. -- 2.39.5