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

a partial template specialization of Functional with test == double, implying "quantity of interest" More...

Public Member Functions

 Functional (std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
 Constructs using a mfem::ParFiniteElementSpace object corresponding to the trial space. More...
 
template<int dim, int... args, typename lambda , typename qpt_data_type = Nothing>
void AddDomainIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, mfem::Mesh &mesh, std::shared_ptr< QuadratureData< qpt_data_type >> qdata=NoQData)
 Adds a domain integral term to the Functional object. 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 , typename qpt_data_type = void>
void AddBoundaryIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, mfem::Mesh &mesh)
 Adds a boundary integral term to the Functional object. 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. 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. 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);
 
double ActionOfGradient (const mfem::Vector &input_T, uint32_t which) const
 this function computes the directional derivative of the quantity of interest functional 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.
 

Detailed Description

template<typename... trials, ExecutionSpace exec>
class serac::Functional< double(trials...), exec >

a partial template specialization of Functional with test == double, implying "quantity of interest"

Definition at line 59 of file functional_qoi.inl.

Constructor & Destructor Documentation

◆ Functional()

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

Constructs using a mfem::ParFiniteElementSpace object corresponding to the trial space.

Parameters
[in]trial_fesThe trial space

Definition at line 83 of file functional_qoi.inl.

Member Function Documentation

◆ ActionOfGradient()

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

this function computes the directional derivative of the quantity of interest functional

Parameters
input_Ta T-vector to apply the action of gradient to
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 282 of file functional_qoi.inl.

◆ AddAreaIntegral()

template<typename... trials, ExecutionSpace exec>
template<int... args, typename lambda , typename qpt_data_type = Nothing>
void serac::Functional< double(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.

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]dataThe data structure containing per-quadrature-point data

Definition at line 242 of file functional_qoi.inl.

◆ AddBoundaryIntegral()

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

Adds a boundary integral term to the Functional object.

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]meshThe 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 204 of file functional_qoi.inl.

◆ AddDomainIntegral()

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

Adds a domain integral term to the Functional object.

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]meshThe domain on which to evaluate the integral
[in]qdataThe data structure containing per-quadrature-point data
Note
The Dimension parameters are used to assist in the deduction of the geometry_dim and spatial_dim template parameter

Definition at line 160 of file functional_qoi.inl.

◆ AddVolumeIntegral()

template<typename... trials, ExecutionSpace exec>
template<int... args, typename lambda , typename qpt_data_type = Nothing>
void serac::Functional< double(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.

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]dataThe data structure containing per-quadrature-point data

Definition at line 259 of file functional_qoi.inl.

◆ operator()()

template<typename... trials, ExecutionSpace exec>
template<uint32_t wrt, typename... T>
operator_paren_return<wrt>::type serac::Functional< double(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

Parameters
tthe time
argsthe input T-vectors

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 323 of file functional_qoi.inl.


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