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

This is a small wrapper around serac::Functional for shape-displaced domains of integration. More...

#include <shape_aware_functional.hpp>

Public Member Functions

template<typename test_type = test, typename = std::enable_if_t<!std::is_same_v<double, test_type>>>
 ShapeAwareFunctional (const mfem::ParFiniteElementSpace *shape_fes, 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<typename test_type = test, typename = std::enable_if_t<std::is_same_v<double, test_type>>>
 ShapeAwareFunctional (const mfem::ParFiniteElementSpace *shape_fes, std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces > trial_fes)
 Constructs a QOI functional using mfem::ParFiniteElementSpace objects corresponding to the trial spaces. More...
 
template<int dim, int... args, typename lambda , typename domain_type , typename qpt_data_type = Nothing>
void AddDomainIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, domain_type &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 domain_type >
void AddBoundaryIntegral (Dimension< dim >, DependsOn< args... >, lambda &&integrand, domain_type &domain)
 Adds a boundary integral term to the weak formulation of the PDE. More...
 
template<uint32_t wrt, typename... T>
auto operator() (DifferentiateWRT< wrt >, double t, const T &... args)
 This function lets the user evaluate the serac::ShapeAwareFunctional with the given trial space values. More...
 
template<typename... T>
auto operator() (double t, const T &... args)
 This function lets the user evaluate the serac::ShapeAwareFunctional with the given trial space values. More...
 
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 shape, typename... trials, ExecutionSpace exec>
class serac::ShapeAwareFunctional< shape, test(trials...), exec >

This is a small wrapper around serac::Functional for shape-displaced domains of integration.

If a finite element kernel is defined on a domain \(x = X + p\) where \(X\) is the reference configuration and \(p\) is a shape displacement field. This shape displacement is typically a input parameter used in shape optimization problems. This wrapper correctly modifies domain and boundary integrals for L2 and H1 trial functions and H1, L2, and QOI (double) test functions for vector-valued H1 shape displacement fields.

See serac::Functional for more details about the underlying finite element abstraction and serac::HeatTransfer for an example of how to use this wrapper.

Template Parameters
testThe space of test function to use
shapeThe space of the shape displacement function to use
trialsThe space of trial functions to use
execwhether to carry out calculations on CPU or GPU

Definition at line 277 of file shape_aware_functional.hpp.

Constructor & Destructor Documentation

◆ ShapeAwareFunctional() [1/2]

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
template<typename test_type = test, typename = std::enable_if_t<!std::is_same_v<double, test_type>>>
serac::ShapeAwareFunctional< shape, test(trials...), exec >::ShapeAwareFunctional ( const mfem::ParFiniteElementSpace *  shape_fes,
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.

Template Parameters
test_typeThe finite element space of the test functions (used to disable this method if called with test = double)
Parameters
[in]shape_fesThe shape displacement finite element space
[in]test_fesThe (non-qoi) test finite element space
[in]trial_fesThe trial finite element spaces

Definition at line 302 of file shape_aware_functional.hpp.

◆ ShapeAwareFunctional() [2/2]

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
template<typename test_type = test, typename = std::enable_if_t<std::is_same_v<double, test_type>>>
serac::ShapeAwareFunctional< shape, test(trials...), exec >::ShapeAwareFunctional ( const mfem::ParFiniteElementSpace *  shape_fes,
std::array< const mfem::ParFiniteElementSpace *, num_trial_spaces >  trial_fes 
)
inline

Constructs a QOI functional using mfem::ParFiniteElementSpace objects corresponding to the trial spaces.

Template Parameters
test_typeThe finite element space of the test functions (used to disable this method if called with test != double)
Parameters
[in]shape_fesThe shape displacement finite element space
[in]trial_fesThe trial finite element spaces

Definition at line 328 of file shape_aware_functional.hpp.

Member Function Documentation

◆ AddBoundaryIntegral()

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

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

Template Parameters
dimThe dimension of the element (2 for quad, 3 for hex, etc)
argsThe type of the trial function input arguments
lambdaThe type of the integrand functor: must implement operator() with an appropriate function signature
domain_typeThe type of the integration domain (either serac::Domain or mfem::Mesh)
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 412 of file shape_aware_functional.hpp.

◆ AddDomainIntegral()

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
template<int dim, int... args, typename lambda , typename domain_type , typename qpt_data_type = Nothing>
void serac::ShapeAwareFunctional< shape, test(trials...), exec >::AddDomainIntegral ( Dimension< dim >  ,
DependsOn< args... >  ,
lambda &&  integrand,
domain_type &  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)
argsThe type of the trial function input arguments
lambdaThe type of the integrand functor: must implement operator() with an appropriate function signature
domain_typeThe type of the integration domain (either serac::Domain or mfem::Mesh)
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
[in,out]qdataThe data for each quadrature point
Note
The Dimension parameters are used to assist in the deduction of the geometry_dim and spatial_dim template parameter

Definition at line 361 of file shape_aware_functional.hpp.

◆ operator()() [1/2]

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
template<uint32_t wrt, typename... T>
auto serac::ShapeAwareFunctional< shape, test(trials...), exec >::operator() ( DifferentiateWRT< wrt >  ,
double  t,
const T &...  args 
)
inline

This function lets the user evaluate the serac::ShapeAwareFunctional with the given trial space values.

Note
The first argument after time in the argument list is always the shape displacement field.
Template Parameters
TThe types of the arguments passed in
wrtThe index of the argument to take the gradient with respect to. Index 0 is always the shape displacement
Parameters
tThe time
argsThe trial space dofs used to carry out the calculation. The first argument is always the shape displacement.
Returns
Either the evaluated integral value or a tuple of the integral value and the requested derivative.

Definition at line 440 of file shape_aware_functional.hpp.

◆ operator()() [2/2]

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
template<typename... T>
auto serac::ShapeAwareFunctional< shape, test(trials...), exec >::operator() ( double  t,
const T &...  args 
)
inline

This function lets the user evaluate the serac::ShapeAwareFunctional with the given trial space values.

Note
The first argument after time in the argument list is always the shape displacement field.
Template Parameters
TThe types of the arguments passed in
Parameters
tThe time
argsThe trial space dofs used to carry out the calculation. The first argument is always the shape displacement. To compute a derivative, at most one argument can be of type differentiate_wrt_this(mfem::Vector).
Returns
Either the evaluated integral value or a tuple of the integral value and the requested derivative.

Definition at line 459 of file shape_aware_functional.hpp.

◆ updateQdata()

template<typename test , typename shape , typename... trials, ExecutionSpace exec>
void serac::ShapeAwareFunctional< shape, 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 472 of file shape_aware_functional.hpp.


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