Reference documentation for deal.II version 9.3.3
\(\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\}}\)
Public Member Functions | Private Attributes | List of all members
DataPostprocessorTensor< dim > Class Template Reference

#include <deal.II/numerics/data_postprocessor.h>

Inheritance diagram for DataPostprocessorTensor< dim >:
[legend]

Public Member Functions

 DataPostprocessorTensor (const std::string &name, const UpdateFlags update_flags)
 
virtual std::vector< std::string > get_names () const override
 
virtual std::vector< DataComponentInterpretation::DataComponentInterpretationget_data_component_interpretation () const override
 
virtual UpdateFlags get_needed_update_flags () const override
 
virtual void evaluate_scalar_field (const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
 
virtual void evaluate_vector_field (const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
 

Private Attributes

const std::string name
 
const UpdateFlags update_flags
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

Detailed Description

template<int dim>
class DataPostprocessorTensor< dim >

This class provides a simpler interface to the functionality offered by the DataPostprocessor class in case one wants to compute only a single tensor quantity (defined as having exactly dim*dim components) from the finite element field passed to the DataOut class.

For this case, we would like to output all of these components as parts of a tensor-valued quantity. Unfortunately, the various backends that write DataOut data in graphical file formats (see the DataOutBase namespace for what formats can be written) do not support tensor data at the current time. In fact, neither does the DataComponentInterpretation namespace that provides semantic information how individual components of graphical data should be interpreted. Nevertheless, like DataPostprocessorScalar and DataPostprocessorVector, this class helps with setting up what the get_names() and get_needed_update_flags() functions required by the DataPostprocessor base class should return, and so the current class implements these based on information that the constructor of the current class receives from further derived classes.

(In order to visualize this collection of scalar fields that, together, are then supposed to be interpreted as a tensor, one has to (i) use a visualization program that can visualize tensors, and (ii) teach it how to re-combine the scalar fields into tensors. In the case of VisIt – see https://wci.llnl.gov/simulation/computer-codes/visit/ – this is done by creating a new "Expression": in essence, one creates a variable, say "grad_u", that is tensor-valued and whose value is given by the expression {{grad_u_xx,grad_u_xy}, {grad_u_yx, grad_u_yy}}, where the referenced variables are the names of scalar fields that, here, are produced by the example below. VisIt is then able to visualize this "new" variable as a tensor.)

All derived classes have to do is implement a constructor and overload either DataPostprocessor::evaluate_scalar_field() or DataPostprocessor::evaluate_vector_field() as discussed in the DataPostprocessor class's documentation.

An example of how the closely related class DataPostprocessorScalar is used can be found in step-29. An example of how the DataPostprocessorVector class can be used is found in the documentation of that class.

An example

A common example of what one wants to do with postprocessors is to visualize not just the value of the solution, but the gradient. This class is meant for tensor-valued outputs, so we will start with a vector-valued solution: the displacement field of step-8. The gradient is a rank-2 tensor (with exactly dim*dim components), so the current class fits the bill to produce the gradient through postprocessing. Then, the following code snippet implements everything you need to have to visualize the gradient:

template <int dim>
class GradientPostprocessor : public DataPostprocessorTensor<dim>
{
public:
GradientPostprocessor ()
:
DataPostprocessorTensor<dim> ("grad_u",
{}
virtual
void
std::vector<Vector<double> > &computed_quantities) const override
{
// ensure that there really are as many output slots
// as there are points at which DataOut provides the
// gradients:
AssertDimension (input_data.solution_gradients.size(),
computed_quantities.size());
for (unsigned int p=0; p<input_data.solution_gradients.size(); ++p)
{
// ensure that each output slot has exactly 'dim*dim'
// components (as should be expected, given that we
// want to create tensor-valued outputs), and copy the
// gradients of the solution at the evaluation points
// into the output slots:
AssertDimension (computed_quantities[p].size(),
for (unsigned int d=0; d<dim; ++d)
for (unsigned int e=0; e<dim; ++e)
= input_data.solution_gradients[p][d][e];
}
}
};
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
Definition: tensor.h:472
Definition: vector.h:110
@ update_gradients
Shape function gradients.
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< std::vector< Tensor< 1, spacedim > > > solution_gradients

The only tricky part in this piece of code is how to sort the dim*dim elements of the strain tensor into the one vector of computed output quantities – in other words, how to unroll the elements of the tensor into the vector. This is facilitated by the Tensor::component_to_unrolled_index() function that takes a pair of indices that specify a particular element of the tensor and returns a vector index that is then used in the code above to fill the computed_quantities array.

The last thing that is necessary is to add another output to the call of DataOut::add_vector() in the output_results() function of the Step8 class of that example program. The corresponding code snippet would then look like this:

GradientPostprocessor<dim> grad_u;
DataOut<dim> data_out;
data_out.attach_dof_handler (dof_handler);
std::vector<DataComponentInterpretation::DataComponentInterpretation>
data_component_interpretation
data_out.add_data_vector (solution,
std::vector<std::string>(dim,"displacement"),
data_component_interpretation);
data_out.add_data_vector (solution, grad_u);
data_out.build_patches ();
data_out.write_vtk (output);
void attach_dof_handler(const DoFHandlerType &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition: data_out.cc:1085
void write_vtk(std::ostream &out) const

This leads to the following output for the displacement field (i.e., the solution) and the gradients (you may want to compare with the solution shown in the results section of step-8; the current data is generated on a uniform mesh for simplicity):

These pictures show an ellipse representing the gradient tensor at, on average, every tenth mesh point. You may want to read through the documentation of the VisIt visualization program (see https://wci.llnl.gov/simulation/computer-codes/visit/) for an interpretation of how exactly tensors are visualizated.

In elasticity, one is often interested not in the gradient of the displacement, but in the "strain", i.e., the symmetrized version of the gradient \(\varepsilon=\frac 12 (\nabla u + \nabla u^T)\). This is easily facilitated with the following minor modification:

template <int dim>
class StrainPostprocessor : public DataPostprocessorTensor<dim>
{
public:
StrainPostprocessor ()
:
DataPostprocessorTensor<dim> ("strain",
{}
virtual
void
std::vector<Vector<double> > &computed_quantities) const override
{
AssertDimension (input_data.solution_gradients.size(),
computed_quantities.size());
for (unsigned int p=0; p<input_data.solution_gradients.size(); ++p)
{
AssertDimension (computed_quantities[p].size(),
for (unsigned int d=0; d<dim; ++d)
for (unsigned int e=0; e<dim; ++e)
= (input_data.solution_gradients[p][d][e]
+
input_data.solution_gradients[p][e][d]) / 2;
}
}
};

Using this class in step-8 leads to the following visualization:

Given how easy it is to output the strain, it would also not be very complicated to write a postprocessor that computes the stress in the solution field as the stress is easily computed from the strain by multiplication with either the strain-stress tensor or, in simple cases, the Lamé constants.

Definition at line 1137 of file data_postprocessor.h.

Constructor & Destructor Documentation

◆ DataPostprocessorTensor()

template<int dim>
DataPostprocessorTensor< dim >::DataPostprocessorTensor ( const std::string &  name,
const UpdateFlags  update_flags 
)

Constructor. Take the name of the single vector variable computed by classes derived from the current one, as well as the update flags necessary to compute this quantity.

Parameters
nameThe name by which the vector variable computed by this class should be made available in graphical output files.
update_flagsThis has to be a combination of update_values, update_gradients, update_hessians and update_quadrature_points. Note that the flag update_quadrature_points updates DataPostprocessorInputs::CommonInputs::evaluation_points. If the DataPostprocessor is to be used in combination with DataOutFaces, you may also ask for a update of normals via the update_normal_vectors flag. The description of the flags can be found at UpdateFlags.

Definition at line 138 of file data_postprocessor.cc.

Member Function Documentation

◆ get_names()

template<int dim>
std::vector< std::string > DataPostprocessorTensor< dim >::get_names
overridevirtual

Return the vector of strings describing the names of the computed quantities. Given the purpose of this class, this is a vector with dim entries all equal to the name given to the constructor.

Implements DataPostprocessor< dim >.

Definition at line 149 of file data_postprocessor.cc.

◆ get_data_component_interpretation()

template<int dim>
std::vector< DataComponentInterpretation::DataComponentInterpretation > DataPostprocessorTensor< dim >::get_data_component_interpretation
overridevirtual

This function returns information about how the individual components of output files that consist of more than one data set are to be interpreted. Since the current class is meant to be used for a single vector result variable, the returned value is obviously DataComponentInterpretation::component_is_part repeated dim times.

Reimplemented from DataPostprocessor< dim >.

Definition at line 158 of file data_postprocessor.cc.

◆ get_needed_update_flags()

template<int dim>
UpdateFlags DataPostprocessorTensor< dim >::get_needed_update_flags
overridevirtual

Return which data has to be provided to compute the derived quantities. The flags returned here are the ones passed to the constructor of this class.

Implements DataPostprocessor< dim >.

Definition at line 167 of file data_postprocessor.cc.

◆ evaluate_scalar_field()

template<int dim>
void DataPostprocessor< dim >::evaluate_scalar_field ( const DataPostprocessorInputs::Scalar< dim > &  input_data,
std::vector< Vector< double > > &  computed_quantities 
) const
virtualinherited

This is the main function which actually performs the postprocessing. The second argument is a reference to the postprocessed data which already has correct size and must be filled by this function.

The function takes the values, gradients, and higher derivatives of the solution at all evaluation points, as well as other data such as the cell, via the first argument. Not all of the member vectors of this argument will be filled with data – in fact, derivatives and other quantities will only be contain valid data if the corresponding flags are returned by an overridden version of the get_needed_update_flags() function (implemented in a user's derived class). Otherwise those vectors will be in an unspecified state.

This function is called when the finite element field that is being converted into graphical data by DataOut or similar classes represents scalar data, i.e., if the finite element in use has only a single real-valued vector component.

Definition at line 26 of file data_postprocessor.cc.

◆ evaluate_vector_field()

template<int dim>
void DataPostprocessor< dim >::evaluate_vector_field ( const DataPostprocessorInputs::Vector< dim > &  input_data,
std::vector< Vector< double > > &  computed_quantities 
) const
virtualinherited

Same as the evaluate_scalar_field() function, but this function is called when the original data vector represents vector data, i.e., the finite element in use has multiple vector components. This function is also called if the finite element is scalar but the solution vector is complex-valued. If the solution vector to be visualized is complex-valued (whether scalar or not), then the input data contains first all real parts of the solution vector at each evaluation point, and then all imaginary parts.

Definition at line 37 of file data_postprocessor.cc.

Member Data Documentation

◆ name

template<int dim>
const std::string DataPostprocessorTensor< dim >::name
private

Copies of the two arguments given to the constructor of this class.

Definition at line 1188 of file data_postprocessor.h.

◆ update_flags

template<int dim>
const UpdateFlags DataPostprocessorTensor< dim >::update_flags
private

Definition at line 1189 of file data_postprocessor.h.


The documentation for this class was generated from the following files: