Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Classes | Public Member Functions | List of all members
serac::Functional< test(trials...), exec > Class Template Reference

Intended to be like std::function for finite element kernels. More...

#include <functional.hpp>

Public Member Functions

 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...
 

Detailed Description

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< H1(Hcurl) > my_residual;
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.
Definition: functional.hpp:209
Template Parameters
testThe space of test functions to use
trialThe space of trial functions to use
execwhether 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
// for domains made up of quadrilaterals embedded in R^2
my_residual.AddAreaIntegral(integrand, domain_of_integration);
// alternatively...
my_residual.AddDomainIntegral(Dimension<2>{}, integrand, domain_of_integration);
// for domains made up of quadrilaterals embedded in R^3
my_residual.AddSurfaceIntegral(integrand, domain_of_integration);
// for domains made up of hexahedra embedded in R^3
my_residual.AddVolumeIntegral(integrand, domain_of_integration);
// alternatively...
my_residual.AddDomainIntegral(Dimension<3>{}, integrand, domain_of_integration);

Definition at line 180 of file functional.hpp.

Constructor & Destructor Documentation

◆ Functional()

template<typename test , typename... trials, ExecutionSpace exec>
serac::Functional< test(trials...), exec >::Functional ( const mfem::ParFiniteElementSpace *  test_fes,
std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces >  trial_fes 
)
inline

Constructs using mfem::ParFiniteElementSpace objects corresponding to the test/trial spaces.

Parameters
[in]test_fesThe (non-qoi) test space
[in]trial_fesThe trial space

Definition at line 209 of file functional.hpp.

Member Function Documentation

◆ ActionOfGradient()

template<typename test , typename... trials, ExecutionSpace exec>
void serac::Functional< test(trials...), exec >::ActionOfGradient ( const mfem::Vector &  input_T,
mfem::Vector &  output_T,
uint32_t  which 
) const
inline

this function computes the directional derivative of serac::Functional::operator()

Parameters
input_Tthe T-vector to apply the action of gradient to
output_Tthe T-vector where the resulting values are stored
whichdescribes which trial space input_T corresponds to
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

Definition at line 394 of file functional.hpp.

◆ AddAreaIntegral()

template<typename test , typename... trials, ExecutionSpace exec>
template<int... args, typename lambda , typename qpt_data_type = Nothing>
void serac::Functional< test(trials...), exec >::AddAreaIntegral ( DependsOn< args... >  which_args,
lambda &&  integrand,
mfem::Mesh &  domain,
std::shared_ptr< QuadratureData< qpt_data_type >>  data = NoQData 
)
inline

Adds an area integral, i.e., over 2D elements in R^2 space.

Template Parameters
lambdathe type of the integrand functor: must implement operator() with an appropriate function signature
qpt_data_typeThe type of the data to store for each quadrature point
Parameters
[in]which_argsa tag type used to indicate which trial spaces are required by this calculation
[in]integrandThe quadrature function
[in]domainThe mesh to evaluate the integral on
[in,out]dataThe data for each quadrature point

Definition at line 354 of file functional.hpp.

◆ AddBoundaryIntegral()

template<typename test , typename... trials, ExecutionSpace exec>
template<int dim, int... args, typename lambda >
void serac::Functional< test(trials...), exec >::AddBoundaryIntegral ( Dimension< dim >  ,
DependsOn< args... >  ,
lambda &&  integrand,
mfem::Mesh &  domain 
)
inline

Adds a boundary integral term to the weak formulation of the PDE.

Template Parameters
dimThe dimension of the boundary element (1 for line, 2 for quad, etc)
lambdathe type of the integrand functor: must implement operator() with an appropriate function signature
Parameters
[in]integrandThe user-provided quadrature function, see Integral
[in]domainThe domain on which to evaluate the integral
Note
The Dimension parameters are used to assist in the deduction of the geometry_dim and spatial_dim template parameter

Definition at line 317 of file functional.hpp.

◆ AddDomainIntegral()

template<typename test , typename... trials, ExecutionSpace exec>
template<int dim, int... args, typename lambda , typename qpt_data_type = Nothing>
void serac::Functional< test(trials...), exec >::AddDomainIntegral ( Dimension< dim >  ,
DependsOn< args... >  ,
lambda &&  integrand,
mfem::Mesh &  domain,
std::shared_ptr< QuadratureData< qpt_data_type >>  qdata = NoQData 
)
inline

Adds a domain integral term to the weak formulation of the PDE.

Template Parameters
dimThe dimension of the element (2 for quad, 3 for hex, etc)
lambdathe type of the integrand functor: must implement operator() with an appropriate function signature
qpt_data_typeThe type of the data to store for each quadrature point
Parameters
[in]integrandThe user-provided quadrature function, see Integral
[in]domainThe domain on which to evaluate the integral
Note
The Dimension parameters are used to assist in the deduction of the geometry_dim and spatial_dim template parameter
Parameters
[in,out]qdataThe data for each quadrature point

Definition at line 275 of file functional.hpp.

◆ AddVolumeIntegral()

template<typename test , typename... trials, ExecutionSpace exec>
template<int... args, typename lambda , typename qpt_data_type = Nothing>
void serac::Functional< test(trials...), exec >::AddVolumeIntegral ( DependsOn< args... >  which_args,
lambda &&  integrand,
mfem::Mesh &  domain,
std::shared_ptr< QuadratureData< qpt_data_type >>  data = NoQData 
)
inline

Adds a volume integral, i.e., over 3D elements in R^3 space.

Template Parameters
lambdathe type of the integrand functor: must implement operator() with an appropriate function signature
qpt_data_typeThe type of the data to store for each quadrature point
Parameters
[in]which_argsa tag type used to indicate which trial spaces are required by this calculation
[in]integrandThe quadrature function
[in]domainThe mesh to evaluate the integral on
[in,out]dataThe data for each quadrature point

Definition at line 370 of file functional.hpp.

◆ operator()()

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
Tthe types of the arguments passed in
Parameters
tthe time
argsthe 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.

◆ updateQdata()

template<typename test , typename... trials, ExecutionSpace exec>
void serac::Functional< test(trials...), exec >::updateQdata ( bool  update_flag)
inline

A flag to update the quadrature data for this operator following the computation.

Typically this is set to false during nonlinear solution iterations and is set to true for the final pass once equilibrium is found.

Parameters
update_flagA flag to update the related quadrature data

Definition at line 511 of file functional.hpp.


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