deal2lkit: A ToolKit library for Deal.II
DOFUtilities Namespace Reference

Define namespace "DOFUtilities" some utilities aimed at solving non-linear problems are implemented. More...

Functions

template<typename Number , typename VEC >
void extract_local_dofs (const VEC &global_vector, const std::vector< types::global_dof_index > &local_dof_indices, std::vector< Number > &independent_local_dofs)
 Extract local dofs values and initialize the number of independent variables up to the second order derivative. More...
 
template<int dim, int spacedim, typename Number >
void get_values (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, std::vector< std::vector< Number > > &us)
 Extract independent local dofs values in a cell and store them in std::vector <std::vector<Number> > us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_values (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Tensor< 1, spacedim, Number > > &us)
 Extract independent local dofs values of a vector_variable in a cell and store them in std::vector <Tensor <1, spacedim, Number> > us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_values (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Scalar &scalar_variable, std::vector< Number > &us)
 Extract independent local dofs values of a scalar_variable in a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_divergences (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Number > &us)
 Extract divergence values of a vector_variable on face of a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_gradients (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Tensor< 2, spacedim, Number > > &grad_us)
 Extract gradient values of a vector_variable on face of a cell and store them in std::vector <Tensor <2, spacedim, Number> > grad_us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_curls (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Tensor< 1,(spacedim > 2 ? spacedim :1), Number > > &curl_us)
 Extract curl values of a vector_variable on face or cell and store them in std::vector <Tensor <1, spacedim, Number> > curl_us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_gradients (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Scalar &scalar_variable, std::vector< Tensor< 1, spacedim, Number > > &grad_us)
 Extract gradient values of a scalar_variable on face of a cell and store them in std::vector <Tensor <1, spacedim, Number> > grad_us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_deformation_gradients (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Tensor< 2, spacedim, Number > > &Fs)
 Compute the deformation gradient in each quadrature point of a cell and it is stored in std::vector <Tensor <2, spacedim, Number> > Fs whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_symmetric_gradients (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &vector_variable, std::vector< Tensor< 2, spacedim, Number > > &grad_us)
 Extract symmetric gradient values of a vector_variable on face of a cell and store them in std::vector <Tensor <2, spacedim, Number> > grad_us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_laplacians (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Scalar &variable, std::vector< Number > &us)
 Extract laplacians local dofs values of a scalar variable in a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_laplacians (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &variable, std::vector< Tensor< 1, spacedim, Number > > &us)
 Extract laplacians local dofs values of a vector variable in a cell and store them in std::vector <Tensor<1,Number > > us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_hessians (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Scalar &variable, std::vector< Tensor< 2, spacedim, Number > > &us)
 Extract hessians local dofs values of a scalar variable in a cell and store them in std::vector <Tensor<2,Number > > us whose size is number of quadrature points in the cell. More...
 
template<int dim, int spacedim, typename Number >
void get_hessians (const FEValuesBase< dim, spacedim > &fe_values, const std::vector< Number > &independent_local_dof_values, const FEValuesExtractors::Vector &variable, std::vector< Tensor< 3, spacedim, Number > > &us)
 Extract hessians local dofs values of a vector variable in a cell and store them in std::vector <Tensor<3,Number > > us whose size is number of quadrature points in the cell. More...
 

Detailed Description

Define namespace "DOFUtilities" some utilities aimed at solving non-linear problems are implemented.

Function Documentation

§ extract_local_dofs()

template<typename Number , typename VEC >
void DOFUtilities::extract_local_dofs ( const VEC &  global_vector,
const std::vector< types::global_dof_index > &  local_dof_indices,
std::vector< Number > &  independent_local_dofs 
)

Extract local dofs values and initialize the number of independent variables up to the second order derivative.

Definition at line 59 of file dof_utilities.h.

§ get_curls()

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_curls ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Tensor< 1,(spacedim > 2 ? spacedim :1), Number >  ,
curl_us 
)

Extract curl values of a vector_variable on face or cell and store them in std::vector <Tensor <1, spacedim, Number> > curl_us whose size is number of quadrature points in the cell.

Definition at line 261 of file dof_utilities.h.

§ get_deformation_gradients()

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_deformation_gradients ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Tensor< 2, spacedim, Number > > &  Fs 
)

Compute the deformation gradient in each quadrature point of a cell and it is stored in std::vector <Tensor <2, spacedim, Number> > Fs whose size is number of quadrature points in the cell.

Definition at line 328 of file dof_utilities.h.

§ get_divergences()

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_divergences ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Number > &  us 
)

Extract divergence values of a vector_variable on face of a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell.

Definition at line 198 of file dof_utilities.h.

§ get_gradients() [1/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_gradients ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Tensor< 2, spacedim, Number > > &  grad_us 
)

Extract gradient values of a vector_variable on face of a cell and store them in std::vector <Tensor <2, spacedim, Number> > grad_us whose size is number of quadrature points in the cell.

Definition at line 228 of file dof_utilities.h.

§ get_gradients() [2/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_gradients ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Scalar scalar_variable,
std::vector< Tensor< 1, spacedim, Number > > &  grad_us 
)

Extract gradient values of a scalar_variable on face of a cell and store them in std::vector <Tensor <1, spacedim, Number> > grad_us whose size is number of quadrature points in the cell.

Definition at line 296 of file dof_utilities.h.

§ get_hessians() [1/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_hessians ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Scalar variable,
std::vector< Tensor< 2, spacedim, Number > > &  us 
)

Extract hessians local dofs values of a scalar variable in a cell and store them in std::vector <Tensor<2,Number > > us whose size is number of quadrature points in the cell.

Definition at line 444 of file dof_utilities.h.

§ get_hessians() [2/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_hessians ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector variable,
std::vector< Tensor< 3, spacedim, Number > > &  us 
)

Extract hessians local dofs values of a vector variable in a cell and store them in std::vector <Tensor<3,Number > > us whose size is number of quadrature points in the cell.

Definition at line 476 of file dof_utilities.h.

§ get_laplacians() [1/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_laplacians ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Scalar variable,
std::vector< Number > &  us 
)

Extract laplacians local dofs values of a scalar variable in a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell.

Definition at line 382 of file dof_utilities.h.

§ get_laplacians() [2/2]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_laplacians ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector variable,
std::vector< Tensor< 1, spacedim, Number > > &  us 
)

Extract laplacians local dofs values of a vector variable in a cell and store them in std::vector <Tensor<1,Number > > us whose size is number of quadrature points in the cell.

Definition at line 413 of file dof_utilities.h.

§ get_symmetric_gradients()

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_symmetric_gradients ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Tensor< 2, spacedim, Number > > &  grad_us 
)

Extract symmetric gradient values of a vector_variable on face of a cell and store them in std::vector <Tensor <2, spacedim, Number> > grad_us whose size is number of quadrature points in the cell.

Definition at line 353 of file dof_utilities.h.

§ get_values() [1/3]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_values ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
std::vector< std::vector< Number > > &  us 
)

Extract independent local dofs values in a cell and store them in std::vector <std::vector<Number> > us whose size is number of quadrature points in the cell.

Definition at line 102 of file dof_utilities.h.

§ get_values() [2/3]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_values ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Vector vector_variable,
std::vector< Tensor< 1, spacedim, Number > > &  us 
)

Extract independent local dofs values of a vector_variable in a cell and store them in std::vector <Tensor <1, spacedim, Number> > us whose size is number of quadrature points in the cell.

Definition at line 137 of file dof_utilities.h.

§ get_values() [3/3]

template<int dim, int spacedim, typename Number >
void DOFUtilities::get_values ( const FEValuesBase< dim, spacedim > &  fe_values,
const std::vector< Number > &  independent_local_dof_values,
const FEValuesExtractors::Scalar scalar_variable,
std::vector< Number > &  us 
)

Extract independent local dofs values of a scalar_variable in a cell and store them in std::vector <Number> us whose size is number of quadrature points in the cell.

Definition at line 168 of file dof_utilities.h.