Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Member Functions | Private Attributes | List of all members
MeshWorker::ScratchData< dim, spacedim > Class Template Reference

#include <deal.II/meshworker/scratch_data.h>

Public Member Functions

 ScratchData (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default)
 
 ScratchData (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const UpdateFlags &neighbor_update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default, const UpdateFlags &neighbor_face_update_flags=update_default)
 
 ScratchData (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default)
 
 ScratchData (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const UpdateFlags &neighbor_update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default, const UpdateFlags &neighbor_face_update_flags=update_default)
 
 ScratchData (const ScratchData< dim, spacedim > &scratch)
 
GeneralDataStorageget_general_data_storage ()
 
const Mapping< dim, spacedim > & get_mapping () const
 
Methods to work on current cell
const FEValues< dim, spacedim > & reinit (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
 
const FEFaceValues< dim, spacedim > & reinit (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_no)
 
const FEFaceValuesBase< dim, spacedim > & reinit (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no)
 
const FEValuesBase< dim, spacedim > & get_current_fe_values () const
 
const std::vector< Point< spacedim > > & get_quadrature_points () const
 
const std::vector< double > & get_JxW_values () const
 
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors () const
 
const std::vector< types::global_dof_index > & get_local_dof_indices () const
 
Methods to work on neighbor cell
const FEValues< dim, spacedim > & reinit_neighbor (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
 
const FEFaceValues< dim, spacedim > & reinit_neighbor (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_no)
 
const FEFaceValuesBase< dim, spacedim > & reinit_neighbor (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no)
 
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values () const
 
const std::vector< double > & get_neighbor_JxW_values () const
 
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors ()
 
const std::vector< types::global_dof_index > & get_neighbor_dof_indices () const
 
Evaluation of finite element fields and their derivatives on the current cell
template<typename VectorType , typename Number = double>
void extract_local_dof_values (const std::string &global_vector_name, const VectorType &input_vector, const Number dummy=Number(0))
 
template<typename Number = double>
const std::vector< Number > & get_local_dof_values (const std::string &global_vector_name, Number dummy=Number(0)) const
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::value_type > & get_values (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::gradient_type > & get_gradients (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::symmetric_gradient_type > & get_symmetric_gradients (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::divergence_type > & get_divergences (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::curl_type > & get_curls (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::hessian_type > & get_hessians (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType< dim, spacedim, Number, Extractor >::third_derivative_type > & get_third_derivatives (const std::string &global_vector_name, const Extractor &variable, const Number dummy=Number(0))
 

Private Member Functions

template<typename Extractor , typename Number = double>
std::string get_unique_name (const std::string &global_vector_name, const Extractor &variable, const std::string &object_type, const unsigned int size, const Number &exemplar_number) const
 
template<typename Number = double>
std::string get_unique_dofs_name (const std::string &global_vector_name, const unsigned int size, const Number &exemplar_number) const
 

Private Attributes

SmartPointer< const Mapping< dim, spacedim > > mapping
 
SmartPointer< const FiniteElement< dim, spacedim > > fe
 
Quadrature< dim > cell_quadrature
 
Quadrature< dim - 1 > face_quadrature
 
UpdateFlags cell_update_flags
 
UpdateFlags neighbor_cell_update_flags
 
UpdateFlags face_update_flags
 
UpdateFlags neighbor_face_update_flags
 
std::unique_ptr< FEValues< dim, spacedim > > fe_values
 
std::unique_ptr< FEFaceValues< dim, spacedim > > fe_face_values
 
std::unique_ptr< FESubfaceValues< dim, spacedim > > fe_subface_values
 
std::unique_ptr< FEValues< dim, spacedim > > neighbor_fe_values
 
std::unique_ptr< FEFaceValues< dim, spacedim > > neighbor_fe_face_values
 
std::unique_ptr< FESubfaceValues< dim, spacedim > > neighbor_fe_subface_values
 
std::vector< types::global_dof_indexlocal_dof_indices
 
std::vector< types::global_dof_indexneighbor_dof_indices
 
GeneralDataStorage user_data_storage
 
GeneralDataStorage internal_data_storage
 
SmartPointer< FEValuesBase< dim, spacedim > > current_fe_values
 
SmartPointer< FEValuesBase< dim, spacedim > > current_neighbor_fe_values
 

Detailed Description

template<int dim, int spacedim = dim>
class MeshWorker::ScratchData< dim, spacedim >

A helper class to simplify the parallel assembly of linear and non-linear problems, and the evaluation of finite element fields.

This class is a drop in ScratchData to use with the WorkStream::run() function and with the MeshWorker::mesh_loop function().

The ScratchData class has three main goals:

The methods in the section "Methods to work on current cell" initialize on demand internal FEValues, FEFaceValues, and FESubfaceValues objects on the current cell, allowing the use of this class as a single substitute for three different objects used to integrate and query finite element values on cells, faces, and subfaces.

Similarly, the methods in the section "Methods to work on neighbor cell" initialize on demand (different) internal FEValues, FEFaceValues, and FESubfaceValues, on neighbor cells, allowing the use of this class also as a single substitute for the additional three objects you would typically need to integrate on the neighbor cell, and on its faces and subfaces (for example, in discontinuous Galerkin methods).

If you need to retrieve values or gradients of finite element solution vectors, on the cell, face, or subface that has just been initialized with one of the functions in the section "Methods to work on current cell", you can use the methods in the section "Evaluation of finite element fields and their derivatives on the current cell".

An example usage for this class is given by the following snippet of code:

ScratchData<dim,spacedim> scratch(fe,
...
for(const auto &cell : dof_handler.active_cell_iterators())
{
const auto &fe_values = scratch.reinit(cell);
scratch.extract_local_dof_values("current_solution",
current_solution);
scratch.extract_local_dof_values("previous_solution",
previous_solution);
const auto &local_solution_values =
scratch.get_local_dof_values("current_solution",
current_solution);
const auto &local_previous_solution_values =
scratch.get_local_dof_values("previous_solution",
previous_solution);
const auto &previous_solution_values_velocity =
scratch.get_values("previous_solution", velocity);
const auto &previous_solution_values_pressure =
scratch.get_values("previous_solution", pressure);
const auto &solution_values_velocity =
scratch.get_values("solution", velocity);
const auto &solution_values_pressure =
scratch.get_values("solution", pressure);
const auto &solution_symmetric_gradients_velocity =
scratch.get_symmetric_gradients("solution", velocity);
// Do something with the above.
}

The order in which you call functions of this class matters: if you call the ScratchData::reinit() function that takes an active cell iterator, then subsequent calls to methods that internally need an FEValuesBase object will use the internal FEValues object initialized with the given cell to perform their calculations. On the other hand, if you have called the ScratchData::reinit() method that also takes a face index, all subsequent calls to methods that need an FEValuesBase object, will use an internally stored FEFaceValues object, initialized with the cell and face index passed to the ScratchData::reinit() function. The same applies for the ScratchData::reinit() method that takes three arguments: the cell, the face index, and the subface index.

The user code should be structured without interleaving work on cells and work on faces.

Consider, for example, the following snippet of code:

ScratchData<dim,spacedim> scratch(...);
for(const auto &cell : dof_handler.active_cell_iterators())
{
const auto &fe_values = scratch.reinit(cell);
const auto &local_dof_values =
scratch.extract_local_dof_values("solution", solution);
// This will contain all values of the temperature on the cell
// quadrature points
const auto &solution_values_cell =
scratch.get_values("solution", temperature);
// Do something with values on the cell
...
// Now start working on the faces
for(unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
{
const auto &fe_face_values = scratch.reinit(cell, f);
// Now we'll have access to the values of the temperature on the faces
const auto &solution_values_face =
scratch.get_values("solution", temperature);
// Notice how the function call is the same, but the result depends
// on what was the last reinit() function that was called. In this
// case, we called reinit(cell, f), triggering internal work on
// faces.
}
// At this point, we would like to go back and work on cells,
// for example querying the values of the gradients:
const auto &solution_gradients =
scratch.get_gradients("solution", temperature);
// This assertion would be triggered in debug mode
AssertDimension(solution_gradients.size(), quadrature_cell.size());
// However, with the above call the content of solution_gradients is the
// gradient of the temperature on the quadrature points of the last face
// visited in the previous loop.
// If you really want to have the values of the gradients on the cell,
// at this point your only option is to call
scratch.reinit(cell);
// (again) before querying for the gradients:
const auto &solution_gradients_cell =
scratch.get_gradients("solution", temperature);
// The call to reinit(cell), however, is an expensive one. You
// should make sure you only call it once per cell by grouping together
// all queries that concern cells in the same place (right after you
// call the reinit(cell) method).
// A similar argument holds for all calls on each face and subface.
}

When using this class, please cite

@article{SartoriGiulianiBardelloni-2018-a,
Author = {Sartori, Alberto and Giuliani, Nicola and
Bardelloni, Mauro and Heltai, Luca},
Journal = {SoftwareX},
Pages = {318--327},
Title = {{deal2lkit: A toolkit library for high performance
programming in deal.II}},
Doi = {10.1016/j.softx.2018.09.004},
Volume = {7},
Year = {2018}}
Author
Luca Heltai, 2019.

Definition at line 232 of file scratch_data.h.

Constructor & Destructor Documentation

◆ ScratchData() [1/5]

template<int dim, int spacedim>
MeshWorker::ScratchData< dim, spacedim >::ScratchData ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags update_flags,
const Quadrature< dim - 1 > &  face_quadrature = Quadrature<dim - 1>(),
const UpdateFlags face_update_flags = update_default 
)

Create an empty ScratchData object. A SmartPointer pointing to mapping and fe is stored internally. Make sure they live longer than this class instance.

The constructor does not initialize any of the internal FEValues objects. These are initialized the first time one of the reinit() functions is called, using the arguments passed here.

Parameters
mappingThe mapping to use in the internal FEValues objects
feThe finite element
quadratureThe cell quadrature
update_flagsUpdateFlags for the current cell FEValues and neighbor cell FEValues
face_quadratureFace quadrature, used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell
face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell

Definition at line 23 of file scratch_data.cc.

◆ ScratchData() [2/5]

template<int dim, int spacedim>
MeshWorker::ScratchData< dim, spacedim >::ScratchData ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags update_flags,
const UpdateFlags neighbor_update_flags,
const Quadrature< dim - 1 > &  face_quadrature = Quadrature<dim - 1>(),
const UpdateFlags face_update_flags = update_default,
const UpdateFlags neighbor_face_update_flags = update_default 
)

Similar to the other constructor, but this one allows to specify different flags for neighbor cells and faces.

Parameters
mappingThe mapping to use in the internal FEValues objects
feThe finite element
quadratureThe cell quadrature
update_flagsUpdateFlags for the current cell FEValues
neighbor_update_flagsUpdateFlags for the neighbor cell FEValues
face_quadratureFace quadrature, used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell
face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for the current cell
neighbor_face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for the neighbor cell

Definition at line 45 of file scratch_data.cc.

◆ ScratchData() [3/5]

template<int dim, int spacedim>
MeshWorker::ScratchData< dim, spacedim >::ScratchData ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags update_flags,
const Quadrature< dim - 1 > &  face_quadrature = Quadrature<dim - 1>(),
const UpdateFlags face_update_flags = update_default 
)

Same as the other constructor, using the default MappingQ1.

Parameters
feThe finite element
quadratureThe cell quadrature
update_flagsUpdateFlags for the current cell FEValues and neighbor cell FEValues
face_quadratureFace quadrature, used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell
face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell

Definition at line 69 of file scratch_data.cc.

◆ ScratchData() [4/5]

template<int dim, int spacedim>
MeshWorker::ScratchData< dim, spacedim >::ScratchData ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags update_flags,
const UpdateFlags neighbor_update_flags,
const Quadrature< dim - 1 > &  face_quadrature = Quadrature<dim - 1>(),
const UpdateFlags face_update_flags = update_default,
const UpdateFlags neighbor_face_update_flags = update_default 
)

Same as the other constructor, using the default MappingQ1.

Parameters
feThe finite element
quadratureThe cell quadrature
update_flagsUpdateFlags for the current cell FEValues
neighbor_update_flagsUpdateFlags for the neighbor cell FEValues
face_quadratureFace quadrature, used for FEFaceValues and FESubfaceValues for both the current cell and the neighbor cell
face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for the current cell
neighbor_face_update_flagsUpdateFlags used for FEFaceValues and FESubfaceValues for the neighbor cell

Definition at line 86 of file scratch_data.cc.

◆ ScratchData() [5/5]

template<int dim, int spacedim>
MeshWorker::ScratchData< dim, spacedim >::ScratchData ( const ScratchData< dim, spacedim > &  scratch)

Deep copy constructor. FEValues objects are not copied.

Definition at line 107 of file scratch_data.cc.

Member Function Documentation

◆ reinit() [1/3]

template<int dim, int spacedim>
const FEValues< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell)

Initialize the internal FEValues with the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

Definition at line 127 of file scratch_data.cc.

◆ reinit() [2/3]

template<int dim, int spacedim>
const FEFaceValues< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell,
const unsigned int  face_no 
)

Initialize the internal FEFaceValues to use the given face_no on the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

Definition at line 144 of file scratch_data.cc.

◆ reinit() [3/3]

template<int dim, int spacedim>
const FEFaceValuesBase< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  subface_no 
)

Initialize the internal FESubfaceValues to use the given subface_no, on face_no, on the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

If subface_no is numbers::invalid_unsigned_int, the reinit() function that takes only the cell and the face_no is called.

Definition at line 162 of file scratch_data.cc.

◆ get_current_fe_values()

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::get_current_fe_values ( ) const

Get the currently initialized FEValues.

This function will return the internal FEValues if the reinit(cell) function was called last. If the reinit(cell, face_no) function was called, then this function returns the internal FEFaceValues, and if the reinit(cell, face_no, subface_no) function was called (with a valid subface_no argument), it returns the internal FESubfaceValues object.

Definition at line 246 of file scratch_data.cc.

◆ get_quadrature_points()

template<int dim, int spacedim>
const std::vector< Point< spacedim > > & MeshWorker::ScratchData< dim, spacedim >::get_quadrature_points ( ) const

Return the quadrature points of the internal FEValues object.

Definition at line 270 of file scratch_data.cc.

◆ get_JxW_values()

template<int dim, int spacedim>
const std::vector< double > & MeshWorker::ScratchData< dim, spacedim >::get_JxW_values ( ) const

Return the JxW values of the internal FEValues object.

Definition at line 279 of file scratch_data.cc.

◆ get_normal_vectors()

template<int dim, int spacedim>
const std::vector< Tensor< 1, spacedim > > & MeshWorker::ScratchData< dim, spacedim >::get_normal_vectors ( ) const

Return the last computed normal vectors.

Definition at line 297 of file scratch_data.cc.

◆ get_local_dof_indices()

template<int dim, int spacedim>
const std::vector< types::global_dof_index > & MeshWorker::ScratchData< dim, spacedim >::get_local_dof_indices ( ) const

Return the local dof indices for the cell passed the last time the reinit() function was called.

Definition at line 315 of file scratch_data.cc.

◆ reinit_neighbor() [1/3]

template<int dim, int spacedim>
const FEValues< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit_neighbor ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell)

Initialize the internal neighbor FEValues to use the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

Definition at line 187 of file scratch_data.cc.

◆ reinit_neighbor() [2/3]

template<int dim, int spacedim>
const FEFaceValues< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit_neighbor ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell,
const unsigned int  face_no 
)

Initialize the internal FEFaceValues to use the given face_no on the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

Definition at line 204 of file scratch_data.cc.

◆ reinit_neighbor() [3/3]

template<int dim, int spacedim>
const FEFaceValuesBase< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::reinit_neighbor ( const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  subface_no 
)

Initialize the internal FESubfaceValues to use the given subface_no, on face_no, on the given cell, and return a reference to it.

After calling this function, get_current_fe_values() will return the same object of this method, as an FEValuesBase reference.

If subface_no is numbers::invalid_unsigned_int, the reinit() function that takes only the cell and the face_no is called.

Definition at line 222 of file scratch_data.cc.

◆ get_current_neighbor_fe_values()

template<int dim, int spacedim>
const FEValuesBase< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::get_current_neighbor_fe_values ( ) const

Get the currently initialized neighbor FEValues.

This function will return the neighbor FEValues if the reinit_neighbor(cell) function was called last. If the reinit_neighbor(cell, face_no) function was called, then this function returns the internal neighbor FEFaceValues, and if the reinit_neighbor(cell, face_no, subface_no) function was called (with a valid subface_no argument), it returns the internal neighbor FESubfaceValues object.

Definition at line 258 of file scratch_data.cc.

◆ get_neighbor_JxW_values()

template<int dim, int spacedim>
const std::vector< double > & MeshWorker::ScratchData< dim, spacedim >::get_neighbor_JxW_values ( ) const

Return the JxW values of the neighbor FEValues object.

Definition at line 288 of file scratch_data.cc.

◆ get_neighbor_normal_vectors()

template<int dim, int spacedim>
const std::vector< Tensor< 1, spacedim > > & MeshWorker::ScratchData< dim, spacedim >::get_neighbor_normal_vectors ( )

Return the last computed normal vectors on the neighbor.

Definition at line 306 of file scratch_data.cc.

◆ get_neighbor_dof_indices()

template<int dim, int spacedim>
const std::vector< types::global_dof_index > & MeshWorker::ScratchData< dim, spacedim >::get_neighbor_dof_indices ( ) const

Return the local dof indices of the neighbor passed the last time the reinit_neighbor() function was called.

Definition at line 324 of file scratch_data.cc.

◆ get_general_data_storage()

template<int dim, int spacedim>
GeneralDataStorage & MeshWorker::ScratchData< dim, spacedim >::get_general_data_storage ( )

Return a GeneralDataStorage object that can be used to store any amount of data, of any type, which is then made accessible by an identifier string.

Definition at line 333 of file scratch_data.cc.

◆ extract_local_dof_values()

template<int dim, int spacedim = dim>
template<typename VectorType , typename Number = double>
void MeshWorker::ScratchData< dim, spacedim >::extract_local_dof_values ( const std::string &  global_vector_name,
const VectorType &  input_vector,
const Number  dummy = Number(0) 
)

Extract the local dof values associated with the internally initialized cell.

Before you call this function, you have to make sure you have previously called one of the reinit() functions.

At every call of this function, a new vector of dof values is generated and stored internally, unless a previous vector with the same name is found. If this is the case, the content of that vector is overwritten.

If you give a unique global_vector_name, then for each cell you are guaranteed to get a unique vector of independent dofs of the same type as the dummy variable. If you use an automatic differentiation number type (like Sacado::Fad::DFad<double>, Sacado::Fad::DFad<Sacado::Fad::DFad<double>>, etc.) this method will also initialize the independent variables internally, allowing you to perform automatic differentiation.

You can access the extracted local dof values by calling the method get_local_dof_values() with the same global_vector_name argument you passed here.

Notice that using this initialization strategy renders the use of this ScratchData object incompatible with the AD helper classes (since they do their own data management). In particular, it is necessary for the user to manage all of the AD data themselves (both before and after this call).

◆ get_local_dof_values()

template<int dim, int spacedim = dim>
template<typename Number = double>
const std::vector<Number>& MeshWorker::ScratchData< dim, spacedim >::get_local_dof_values ( const std::string &  global_vector_name,
Number  dummy = Number(0) 
) const

After calling extract_local_dof_values(), you can retrieve the stored information through this method.

Both the argument global_vector_name and the type of the dummy variable should match the ones you passed to the extract_local_dof_values() function.

◆ get_values()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: value_type>& MeshWorker::ScratchData< dim, spacedim >::get_values ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the values of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_gradients()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: gradient_type>& MeshWorker::ScratchData< dim, spacedim >::get_gradients ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the gradients of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_symmetric_gradients()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: symmetric_gradient_type>& MeshWorker::ScratchData< dim, spacedim >::get_symmetric_gradients ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the symmetric gradients of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_divergences()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: divergence_type>& MeshWorker::ScratchData< dim, spacedim >::get_divergences ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the divergences of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_curls()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector<typename internal:: OutputType<dim, spacedim, Number, Extractor>::curl_type>& MeshWorker::ScratchData< dim, spacedim >::get_curls ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the curls of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_hessians()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: hessian_type>& MeshWorker::ScratchData< dim, spacedim >::get_hessians ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the hessians of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues class for more information.

◆ get_third_derivatives()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
const std::vector< typename internal::OutputType<dim, spacedim, Number, Extractor>:: third_derivative_type>& MeshWorker::ScratchData< dim, spacedim >::get_third_derivatives ( const std::string &  global_vector_name,
const Extractor &  variable,
const Number  dummy = Number(0) 
)

For the solution vector identified by global_vector_name, compute the third_derivatives of the function at the quadrature points, and return a vector with the correct type deduced by the Extractor you passed as the variable argument.

Before you can call this method, you need to call the extract_local_dof_values() method at least once, passing the same global_vector_name string, and the same type for the variable dummy.

If you have not previously called the extract_local_dof_values() method, this function will throw an exception.

For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object for which you called one of the reinit() functions, must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See "The interplay of UpdateFlags, Mapping, and FiniteElement" in the documentation of the FEValues for more information.

◆ get_mapping()

template<int dim, int spacedim>
const Mapping< dim, spacedim > & MeshWorker::ScratchData< dim, spacedim >::get_mapping ( ) const

Return a reference to the used mapping.

Definition at line 342 of file scratch_data.cc.

◆ get_unique_name()

template<int dim, int spacedim = dim>
template<typename Extractor , typename Number = double>
std::string MeshWorker::ScratchData< dim, spacedim >::get_unique_name ( const std::string &  global_vector_name,
const Extractor &  variable,
const std::string &  object_type,
const unsigned int  size,
const Number &  exemplar_number 
) const
private

Construct a unique name to store vectors of values, gradients, divergences, etc., in the internal GeneralDataStorage object.

◆ get_unique_dofs_name()

template<int dim, int spacedim = dim>
template<typename Number = double>
std::string MeshWorker::ScratchData< dim, spacedim >::get_unique_dofs_name ( const std::string &  global_vector_name,
const unsigned int  size,
const Number &  exemplar_number 
) const
private

Construct a unique name to store local dof values.

Member Data Documentation

◆ mapping

template<int dim, int spacedim = dim>
SmartPointer<const Mapping<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::mapping
private

The mapping used by the internal FEValues. Make sure it lives longer than this class.

Definition at line 796 of file scratch_data.h.

◆ fe

template<int dim, int spacedim = dim>
SmartPointer<const FiniteElement<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::fe
private

The finite element used by the internal FEValues. Make sure it lives longer than this class.

Definition at line 802 of file scratch_data.h.

◆ cell_quadrature

template<int dim, int spacedim = dim>
Quadrature<dim> MeshWorker::ScratchData< dim, spacedim >::cell_quadrature
private

Quadrature formula used to integrate on the current cell, and on its neighbor.

Definition at line 808 of file scratch_data.h.

◆ face_quadrature

template<int dim, int spacedim = dim>
Quadrature<dim - 1> MeshWorker::ScratchData< dim, spacedim >::face_quadrature
private

Quadrature formula used to integrate on faces, subfaces, and neighbor faces and subfaces.

Definition at line 814 of file scratch_data.h.

◆ cell_update_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::ScratchData< dim, spacedim >::cell_update_flags
private

UpdateFlags to use when initializing the cell FEValues object.

Definition at line 819 of file scratch_data.h.

◆ neighbor_cell_update_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::ScratchData< dim, spacedim >::neighbor_cell_update_flags
private

UpdateFlags to use when initializing the neighbor cell FEValues objects.

Definition at line 824 of file scratch_data.h.

◆ face_update_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::ScratchData< dim, spacedim >::face_update_flags
private

UpdateFlags to use when initializing FEFaceValues and FESubfaceValues objects.

Definition at line 830 of file scratch_data.h.

◆ neighbor_face_update_flags

template<int dim, int spacedim = dim>
UpdateFlags MeshWorker::ScratchData< dim, spacedim >::neighbor_face_update_flags
private

UpdateFlags to use when initializing neighbor FEFaceValues and FESubfaceValues objects.

Definition at line 836 of file scratch_data.h.

◆ fe_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::fe_values
private

Finite element values on the current cell.

Definition at line 841 of file scratch_data.h.

◆ fe_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEFaceValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::fe_face_values
private

Finite element values on the current face.

Definition at line 846 of file scratch_data.h.

◆ fe_subface_values

template<int dim, int spacedim = dim>
std::unique_ptr<FESubfaceValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::fe_subface_values
private

Finite element values on the current subface.

Definition at line 851 of file scratch_data.h.

◆ neighbor_fe_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::neighbor_fe_values
private

Finite element values on the neighbor cell.

Definition at line 856 of file scratch_data.h.

◆ neighbor_fe_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEFaceValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::neighbor_fe_face_values
private

Finite element values on the neighbor face.

Definition at line 861 of file scratch_data.h.

◆ neighbor_fe_subface_values

template<int dim, int spacedim = dim>
std::unique_ptr<FESubfaceValues<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::neighbor_fe_subface_values
private

Finite element values on the neighbor subface.

Definition at line 866 of file scratch_data.h.

◆ local_dof_indices

template<int dim, int spacedim = dim>
std::vector<types::global_dof_index> MeshWorker::ScratchData< dim, spacedim >::local_dof_indices
private

Dof indices on the current cell.

Definition at line 871 of file scratch_data.h.

◆ neighbor_dof_indices

template<int dim, int spacedim = dim>
std::vector<types::global_dof_index> MeshWorker::ScratchData< dim, spacedim >::neighbor_dof_indices
private

Dof indices on the neighbor cell.

Definition at line 876 of file scratch_data.h.

◆ user_data_storage

template<int dim, int spacedim = dim>
GeneralDataStorage MeshWorker::ScratchData< dim, spacedim >::user_data_storage
private

User data storage.

Definition at line 881 of file scratch_data.h.

◆ internal_data_storage

template<int dim, int spacedim = dim>
GeneralDataStorage MeshWorker::ScratchData< dim, spacedim >::internal_data_storage
private

Internal data storage.

Definition at line 886 of file scratch_data.h.

◆ current_fe_values

template<int dim, int spacedim = dim>
SmartPointer<FEValuesBase<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::current_fe_values
private

A pointer to the last used FEValues/FEFaceValues, or FESubfaceValues object on this cell.

Definition at line 892 of file scratch_data.h.

◆ current_neighbor_fe_values

template<int dim, int spacedim = dim>
SmartPointer<FEValuesBase<dim, spacedim> > MeshWorker::ScratchData< dim, spacedim >::current_neighbor_fe_values
private

A pointer to the last used FEValues/FEFaceValues, or FESubfaceValues object on the neighbor cell.

Definition at line 898 of file scratch_data.h.


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