a class for representing a Integral calculations and their derivatives
More...
#include <integral.hpp>
|
|
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
|
| |
|
| | 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 &dinput_E, mfem::BlockVector &doutput_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...
|
| |
| bool | DependsOn (uint32_t which) const |
| | returns whether or not this integral depends on argument which More...
|
| |
|
|
static constexpr std::size_t | num_types = 3 |
| | the number of different kinds of integration domains
|
| |
a class for representing a Integral calculations and their derivatives
Definition at line 25 of file integral.hpp.
◆ 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
-
| d | the domain of integration |
| trial_space_indices | a list of which trial spaces are used in the integrand |
Definition at line 38 of file integral.hpp.
◆ 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_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 114 of file integral.hpp.
◆ DependsOn()
| bool serac::Integral::DependsOn |
( |
uint32_t |
which | ) |
const |
|
inline |
returns whether or not this integral depends on argument which
- Parameters
-
- Returns
- true when if this Integral object was created with a DependsOn<...> statement that includes
i
Definition at line 131 of file integral.hpp.
◆ GradientMult()
| void serac::Integral::GradientMult |
( |
const mfem::BlockVector & |
dinput_E, |
|
|
mfem::BlockVector & |
doutput_E, |
|
|
uint32_t |
differentiation_index |
|
) |
| const |
|
inline |
evaluate the jacobian(with respect to some trial space)-vector product of this integral
- Parameters
-
| dinput_E | a block vector (block index corresponds to the element geometry) of a specific trial space element values |
| doutput_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 94 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
-
| 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 66 of file integral.hpp.
◆ 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 180 of file integral.hpp.
The documentation for this struct was generated from the following file: