Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
serac::Integral Struct Reference

a class for representing a Integral calculations and their derivatives More...

#include <integral.hpp>

Collaboration diagram for serac::Integral:
Collaboration graph
[legend]

Public Types

using eval_func = std::function< void(double, const std::vector< const double * > &, double *, bool)>
 signature of integral evaluation kernel
 
using jacobian_vector_product_func = std::function< void(const double *, double *)>
 signature of element jvp kernel
 
using grad_func = std::function< void(ExecArrayView< double, 3, ExecutionSpace::CPU >)>
 signature of element gradient kernel
 

Public Member Functions

 Integral (const Domain &d, std::vector< uint32_t > trial_space_indices)
 Construct an "empty" Integral object, whose kernels are to be initialized later in one of the Make****Integral free functions below. More...
 
void Mult (double t, const std::vector< mfem::BlockVector > &input_E, mfem::BlockVector &output_E, uint32_t differentiation_index, bool update_state) const
 evaluate the integral, optionally storing q-function derivatives with respect to a specific trial space. More...
 
void GradientMult (const mfem::BlockVector &input_E, mfem::BlockVector &output_E, uint32_t differentiation_index) const
 evaluate the jacobian(with respect to some trial space)-vector product of this integral More...
 
void ComputeElementGradients (std::map< mfem::Geometry::Type, ExecArray< double, 3, ExecutionSpace::CPU > > &K_e, uint32_t differentiation_index) const
 evaluate the jacobian (with respect to some trial space) of this integral More...
 

Public Attributes

Domain domain_
 information about which elements to integrate over
 
std::map< mfem::Geometry::Type, eval_funcevaluation_
 kernels for integral evaluation over each type of element
 
std::vector< std::map< mfem::Geometry::Type, eval_func > > evaluation_with_AD_
 kernels for integral evaluation + derivative w.r.t. specified argument over each type of element
 
std::vector< std::map< mfem::Geometry::Type, jacobian_vector_product_func > > jvp_
 kernels for jacobian-vector product of integral calculation
 
std::vector< std::map< mfem::Geometry::Type, grad_func > > element_gradient_
 kernels for calculation of element jacobians
 
std::vector< uint32_t > active_trial_spaces_
 a list of the trial spaces that take part in this integrand
 
std::map< uint32_t, uint32_t > functional_to_integral_index_
 a way of translating between the indices used by Functional and Integral to refer to the same trial space. More...
 
std::map< mfem::Geometry::Type, GeometricFactorsgeometric_factors_
 the spatial positions and jacobians (dx_dxi) for each element type and quadrature point
 

Static Public Attributes

static constexpr std::size_t num_types = 2
 the number of different kinds of integration domains
 

Detailed Description

a class for representing a Integral calculations and their derivatives

Definition at line 24 of file integral.hpp.

Constructor & Destructor Documentation

◆ Integral()

serac::Integral::Integral ( const Domain d,
std::vector< uint32_t >  trial_space_indices 
)
inline

Construct an "empty" Integral object, whose kernels are to be initialized later in one of the Make****Integral free functions below.

Note
It is not intended that users construct these objects manually
Parameters
dthe domain of integration
trial_space_indicesa list of which trial spaces are used in the integrand

Definition at line 37 of file integral.hpp.

Member Function Documentation

◆ ComputeElementGradients()

void serac::Integral::ComputeElementGradients ( std::map< mfem::Geometry::Type, ExecArray< double, 3, ExecutionSpace::CPU > > &  K_e,
uint32_t  differentiation_index 
) const
inline

evaluate the jacobian (with respect to some trial space) of this integral

Parameters
K_ea collection (one for each element type) of element jacobians (num_elements x trial_dofs_per_elem x test_dofs_per_elem)
differentiation_indexthe index of the trial space being differentiated

Definition at line 112 of file integral.hpp.

◆ GradientMult()

void serac::Integral::GradientMult ( const mfem::BlockVector &  input_E,
mfem::BlockVector &  output_E,
uint32_t  differentiation_index 
) const
inline

evaluate the jacobian(with respect to some trial space)-vector product of this integral

Parameters
input_Ea block vector (block index corresponds to the element geometry) of a specific trial space element values
output_Ea block vector (block index corresponds to the element geometry) of the output values for each element.
differentiation_indexa non-negative value indicates directional derivative with respect to the trial space with that index.

Definition at line 93 of file integral.hpp.

◆ Mult()

void serac::Integral::Mult ( double  t,
const std::vector< mfem::BlockVector > &  input_E,
mfem::BlockVector &  output_E,
uint32_t  differentiation_index,
bool  update_state 
) const
inline

evaluate the integral, optionally storing q-function derivatives with respect to a specific trial space.

Parameters
tthe time
input_Ea collection (one for each trial space) of block vectors (block index corresponds to the element geometry) containing input values for each element.
output_Ea block vector (block index corresponds to the element geometry) of the output values for each element.
differentiation_indexa non-negative value indicates differentiation with respect to the trial space with that index. A value of -1 indicates no differentiation will occur.
update_statewhether or not to store the updated state values computed in the q-function. For plasticity and other path-dependent materials, this flag should only be set to true once a solution to the nonlinear system has been found.

Definition at line 65 of file integral.hpp.

Member Data Documentation

◆ functional_to_integral_index_

std::map<uint32_t, uint32_t> serac::Integral::functional_to_integral_index_

a way of translating between the indices used by Functional and Integral to refer to the same trial space.

e.g. A Functional may have 4 trial spaces {A, B, C, D}, but an Integral may only depend on a subset, say {B, C}. From Functional's perspective, trial spaces {B, C} have (zero-based) indices of {1, 2}, but from Integral's perspective, those are trial spaces {0, 1}.

So, in this example functional_to_integral_index_ would have the values:

std::map<int,int> functional_to_integral = {{1, 0}, {2, 1}};

Definition at line 164 of file integral.hpp.


The documentation for this struct was generated from the following file: