|
| Functional (const mfem::ParFiniteElementSpace *test_fes, std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes) |
| Constructs using mfem::ParFiniteElementSpace objects corresponding to the test/trial spaces. More...
|
|
template<int dim, int... args, typename lambda , typename qpt_data_type = Nothing> |
void | AddDomainIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, mfem::Mesh &domain, std::shared_ptr< QuadratureData< qpt_data_type >> qdata=NoQData) |
| Adds a domain integral term to the weak formulation of the PDE. More...
|
|
template<int dim, int... args, typename lambda , typename qpt_data_type = Nothing> |
void | AddDomainIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, Domain &domain, std::shared_ptr< QuadratureData< qpt_data_type >> qdata=NoQData) |
| This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
|
template<int dim, int... args, typename lambda > |
void | AddBoundaryIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, mfem::Mesh &domain) |
| Adds a boundary integral term to the weak formulation of the PDE. More...
|
|
template<int dim, int... args, typename lambda > |
void | AddBoundaryIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, const Domain &domain) |
| This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
|
template<int... args, typename lambda , typename qpt_data_type = Nothing> |
void | AddAreaIntegral (DependsOn< args... > which_args, lambda &&integrand, mfem::Mesh &domain, std::shared_ptr< QuadratureData< qpt_data_type >> data=NoQData) |
| Adds an area integral, i.e., over 2D elements in R^2 space. More...
|
|
template<int... args, typename lambda , typename qpt_data_type = Nothing> |
void | AddVolumeIntegral (DependsOn< args... > which_args, lambda &&integrand, mfem::Mesh &domain, std::shared_ptr< QuadratureData< qpt_data_type >> data=NoQData) |
| Adds a volume integral, i.e., over 3D elements in R^3 space. More...
|
|
template<int... args, typename lambda > |
void | AddSurfaceIntegral (DependsOn< args... > which_args, lambda &&integrand, mfem::Mesh &domain) |
| alias for Functional::AddBoundaryIntegral(Dimension<2>{}, integrand, domain);
|
|
void | ActionOfGradient (const mfem::Vector &input_T, mfem::Vector &output_T, uint32_t which) const |
| this function computes the directional derivative of serac::Functional::operator() More...
|
|
template<uint32_t wrt, typename... T> |
operator_paren_return< wrt >::type | operator() (DifferentiateWRT< wrt >, double t, const T &... args) |
| this function lets the user evaluate the serac::Functional with the given trial space values More...
|
|
template<typename... T> |
auto | operator() (double t, const T &... args) |
| This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
|
void | updateQdata (bool update_flag) |
| A flag to update the quadrature data for this operator following the computation. More...
|
|
template<typename test, typename... trials, ExecutionSpace exec>
class serac::Functional< test(trials...), exec >
Intended to be like std::function
for finite element kernels.
That is: you tell it the inputs (trial spaces) for a kernel, and the outputs (test space) like std::function
.
For example, this code represents a function that takes an integer argument and returns a double:
std::function< double(int) > my_func;
And this represents a function that takes values from an Hcurl field and returns a residual vector associated with an H1 field:
Functional(const mfem::ParFiniteElementSpace *test_fes, std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
Constructs using mfem::ParFiniteElementSpace objects corresponding to the test/trial spaces.
- Template Parameters
-
test | The space of test functions to use |
trial | The space of trial functions to use |
exec | whether to carry out calculations on CPU or GPU |
To use this class, you use the methods Functional::Add****Integral
(integrand,domain_of_integration) where integrand
is a q-function lambda or functor and domain_of_integration
is an mfem::mesh
- See also
- https://libceed.readthedocs.io/en/latest/libCEEDapi/#theoretical-framework for additional information on the idea behind a quadrature function and its inputs/outputs
my_residual.AddAreaIntegral(integrand, domain_of_integration);
my_residual.AddDomainIntegral(Dimension<2>{}, integrand, domain_of_integration);
my_residual.AddSurfaceIntegral(integrand, domain_of_integration);
my_residual.AddVolumeIntegral(integrand, domain_of_integration);
my_residual.AddDomainIntegral(Dimension<3>{}, integrand, domain_of_integration);
Definition at line 180 of file functional.hpp.
template<typename test , typename... trials, ExecutionSpace exec>
template<uint32_t wrt, typename... T>
operator_paren_return<wrt>::type serac::Functional< test(trials...), exec >::operator() |
( |
DifferentiateWRT< wrt > |
, |
|
|
double |
t, |
|
|
const T &... |
args |
|
) |
| |
|
inline |
this function lets the user evaluate the serac::Functional with the given trial space values
note: it accepts exactly num_trial_spaces
arguments of type mfem::Vector. Additionally, one of those arguments may be a dual_vector, to indicate that Functional::operator() should not only evaluate the element calculations, but also differentiate them w.r.t. the specified dual_vector argument
- Template Parameters
-
T | the types of the arguments passed in |
- Parameters
-
t | the time |
args | the trial space dofs used to carry out the calculation, at most one of which may be of the type differentiate_wrt_this(mfem::Vector) |
Definition at line 435 of file functional.hpp.