Serac
0.1
Serac is an implicit thermal strucural mechanics simulation code.
|
a class for representing a Integral calculations and their derivatives More...
#include <integral.hpp>
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_func > | evaluation_ |
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, GeometricFactors > | geometric_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 | |
a class for representing a Integral calculations and their derivatives
Definition at line 24 of file integral.hpp.
|
inline |
Construct an "empty" Integral object, whose kernels are to be initialized later in one of the Make****Integral free functions below.
d | the domain of integration |
trial_space_indices | a list of which trial spaces are used in the integrand |
Definition at line 37 of file integral.hpp.
|
inline |
evaluate the jacobian (with respect to some trial space) of this integral
K_e | a collection (one for each element type) of element jacobians (num_elements x trial_dofs_per_elem x test_dofs_per_elem) |
differentiation_index | the index of the trial space being differentiated |
Definition at line 112 of file integral.hpp.
|
inline |
evaluate the jacobian(with respect to some trial space)-vector product of this integral
input_E | a block vector (block index corresponds to the element geometry) of a specific trial space element values |
output_E | a block vector (block index corresponds to the element geometry) of the output values for each element. |
differentiation_index | a non-negative value indicates directional derivative with respect to the trial space with that index. |
Definition at line 93 of file integral.hpp.
|
inline |
evaluate the integral, optionally storing q-function derivatives with respect to a specific trial space.
t | the time |
input_E | a collection (one for each trial space) of block vectors (block index corresponds to the element geometry) containing input values for each element. |
output_E | a block vector (block index corresponds to the element geometry) of the output values for each element. |
differentiation_index | a non-negative value indicates differentiation with respect to the trial space with that index. A value of -1 indicates no differentiation will occur. |
update_state | whether 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.
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:
Definition at line 164 of file integral.hpp.