Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Classes | Public Member Functions | Static Public Attributes | Protected Types | Protected Attributes | List of all members
serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > > Class Template Reference

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. More...

#include <heat_transfer.hpp>

Collaboration diagram for serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >:
Collaboration graph
[legend]

Classes

struct  ThermalMaterialIntegrand
 Functor representing the integrand of a thermal material. Material type must be a functor as well. More...
 

Public Member Functions

 HeatTransfer (const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::string mesh_tag, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false)
 Construct a new heat transfer object. More...
 
 HeatTransfer (std::unique_ptr< serac::EquationSolver > solver, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::string mesh_tag, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false)
 Construct a new heat transfer object. More...
 
 HeatTransfer (const HeatTransferInputOptions &input_options, const std::string &physics_name, const std::string &mesh_tag, int cycle=0, double time=0.0)
 Construct a new Nonlinear HeatTransfer Solver object. More...
 
void initializeThermalStates ()
 Non virtual method to reset thermal states to zero. This does not reset design parameters or shape. More...
 
void resetStates (int cycle=0, double time=0.0) override
 Method to reset physics states to zero. This does not reset design parameters or shape. More...
 
void setTemperatureBCs (const std::set< int > &temp_bdr, std::function< double(const mfem::Vector &x, double t)> temp)
 Set essential temperature boundary conditions (strongly enforced) More...
 
void advanceTimestep (double dt) override
 Advance the heat conduction physics module in time. More...
 
template<int... active_parameters, typename MaterialType >
void setMaterial (DependsOn< active_parameters... >, MaterialType material)
 Set the thermal material model for the physics solver. More...
 
template<typename MaterialType >
void setMaterial (MaterialType material)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void setTemperature (std::function< double(const mfem::Vector &x, double t)> temp)
 Set the underlying finite element state to a prescribed temperature. More...
 
void setTemperature (const FiniteElementState temp)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<int... active_parameters, typename SourceType >
void setSource (DependsOn< active_parameters... >, SourceType source_function, const std::optional< Domain > &optional_domain=std::nullopt)
 Set the thermal source function. More...
 
template<typename SourceType >
void setSource (SourceType source_function, const std::optional< Domain > &optional_domain=std::nullopt)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<int... active_parameters, typename FluxType >
void setFluxBCs (DependsOn< active_parameters... >, FluxType flux_function, const std::optional< Domain > &optional_domain=std::nullopt)
 Set the thermal flux boundary condition. More...
 
template<typename FluxType >
void setFluxBCs (FluxType flux_function, const std::optional< Domain > &optional_domain=std::nullopt)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void setShapeDisplacement (std::function< void(const mfem::Vector &x, mfem::Vector &shape_disp)> shape_disp)
 Set the underlying finite element state to a prescribed shape displacement. More...
 
void setShapeDisplacement (FiniteElementState &shape_disp)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
const serac::FiniteElementStatetemperature () const
 Get the temperature state. More...
 
const serac::FiniteElementStatetemperatureRate () const
 Get the temperature rate of change state. More...
 
const FiniteElementStatestate (const std::string &state_name) const override
 Accessor for getting named finite element state primal solution from the physics modules. More...
 
void setState (const std::string &state_name, const FiniteElementState &state) override
 Set the primal solution field (temperature) for the underlying heat transfer solver. More...
 
virtual std::vector< std::string > stateNames () const override
 Get a vector of the finite element state primal solution names. More...
 
template<int... active_parameters, typename callable >
void addCustomBoundaryIntegral (DependsOn< active_parameters... >, callable qfunction, const std::optional< Domain > &optional_domain=std::nullopt)
 register a custom boundary integral calculation as part of the residual More...
 
const FiniteElementStateadjoint (const std::string &state_name) const override
 Accessor for getting named finite element state adjoint solution from the physics modules. More...
 
void completeSetup () override
 Complete the initialization and allocation of the data structures. More...
 
virtual void setAdjointLoad (std::unordered_map< std::string, const serac::FiniteElementDual & > loads) override
 Set the loads for the adjoint reverse timestep solve. More...
 
void reverseAdjointTimestep () override
 Solve the adjoint problem. More...
 
std::unordered_map< std::string, FiniteElementStategetCheckpointedStates (int cycle_to_load) const override
 Accessor for getting named finite element state primal solution from the physics modules at a given checkpointed cycle index. More...
 
FiniteElementDualcomputeTimestepSensitivity (size_t parameter_field) override
 Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoint problem with respect to the parameter field for the last computed adjoint timestep. More...
 
FiniteElementDualcomputeTimestepShapeSensitivity () override
 Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoint problem with respect to the shape displacement field for the last computed adjoint timestep. More...
 
const std::unordered_map< std::string, const serac::FiniteElementDual & > computeInitialConditionSensitivity () override
 Compute the implicit sensitivity of the quantity of interest with respect to the initial temperature. More...
 
virtual ~HeatTransfer ()=default
 Destroy the Thermal Solver object.
 
- Public Member Functions inherited from serac::BasePhysics
 BasePhysics (std::string physics_name, std::string mesh_tag, int cycle=0, double time=0.0, bool checkpoint_to_disk=false)
 Empty constructor. More...
 
 BasePhysics (BasePhysics &&other)=default
 Construct a new Base Physics object (copy constructor) More...
 
virtual double time () const
 Get the current forward-solution time. More...
 
virtual int cycle () const
 Get the current forward-solution cycle iteration number. More...
 
virtual double maxTime () const
 Get the maximum time reached by the forward solver. More...
 
virtual double minTime () const
 Get the initial time used by the forward solver. More...
 
virtual int maxCycle () const
 The maximum cycle (timestep iteration number) reached by the forward solver. More...
 
virtual int minCycle () const
 Get the initial cycle (timestep iteration number) used by the forward solver. More...
 
bool isQuasistatic () const
 Check if the physics is setup as quasistatic. More...
 
virtual std::vector< double > timesteps () const
 Get a vector of the timestep sizes (i.e. \(\Delta t\)s) taken by the forward solver. More...
 
virtual std::vector< std::string > adjointNames () const
 Get a vector of the finite element state adjoint solution names. More...
 
const FiniteElementStateshapeDisplacement () const
 Accessor for getting the shape displacement field from the physics modules. More...
 
const FiniteElementStateparameter (const std::string &parameter_name) const
 Accessor for getting named finite element state parameter fields from the physics modules. More...
 
const FiniteElementStateparameter (std::size_t parameter_index) const
 Accessor for getting indexed finite element state parameter fields from the physics modules. More...
 
std::vector< std::string > parameterNames ()
 Get a vector of the finite element state parameter names. More...
 
void setParameter (const size_t parameter_index, const FiniteElementState &parameter_state)
 Deep copy a parameter field into the internally-owned parameter used for simulations. More...
 
void setShapeDisplacement (const FiniteElementState &shape_displacement)
 Set the current shape displacement for the underlying mesh. More...
 
virtual void outputStateToDisk (std::optional< std::string > paraview_output_dir={}) const
 Output the current state of the PDE fields in Sidre format and optionally in Paraview format if paraview_output_dir is given. More...
 
FiniteElementState loadCheckpointedState (const std::string &state_name, int cycle) const
 Accessor for getting a single named finite element state primal solution from the physics modules at a given checkpointed cycle index. More...
 
virtual double getCheckpointedTimestep (int cycle) const
 Get a timestep increment which has been previously checkpointed at the give cycle. More...
 
virtual void initializeSummary (axom::sidre::DataStore &datastore, const double t_final, const double dt) const
 Initializes the Sidre structure for simulation summary data. More...
 
virtual void saveSummary (axom::sidre::DataStore &datastore, const double t) const
 Saves the summary data to the Sidre Datastore. More...
 
virtual ~BasePhysics ()=default
 Destroy the Base Solver object.
 
const mfem::ParMesh & mesh () const
 Returns a reference to the mesh object.
 
mfem::ParMesh & mesh ()
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 

Static Public Attributes

static constexpr auto NUM_STATE_VARS = 2
 The total number of non-parameter state variables (temperature, dtemp_dt) passed to the FEM integrators.
 

Protected Types

using scalar_trial = H1< order >
 The compile-time finite element trial space for heat transfer (H1 of order p)
 
using shape_trial = H1< SHAPE_ORDER, dim >
 The compile-time finite element trial space for shape displacement (vector H1 of order 1)
 
using test = H1< order >
 The compile-time finite element test space for heat transfer (H1 of order p)
 

Protected Attributes

serac::FiniteElementState temperature_
 The temperature finite element state.
 
FiniteElementState temperature_rate_
 Rate of change in temperature at the current adjoint timestep.
 
serac::FiniteElementState adjoint_temperature_
 The adjoint temperature finite element states, the multiplier on the residual for a given timestep.
 
serac::FiniteElementDual implicit_sensitivity_temperature_start_of_step_
 The total/implicit sensitivity of the qoi with respect to the start of the previous timestep's temperature.
 
serac::FiniteElementDual temperature_adjoint_load_
 The downstream derivative of the qoi with repect to the primal temperature variable.
 
serac::FiniteElementDual temperature_rate_adjoint_load_
 The downstream derivative of the qoi with repect to the primal temperature rate variable.
 
std::unique_ptr< ShapeAwareFunctional< shape_trial, test(scalar_trial, scalar_trial, parameter_space...)> > residual_
 serac::Functional that is used to calculate the residual and its derivatives
 
std::unique_ptr< mfem::HypreParMatrix > M_
 Assembled mass matrix.
 
std::shared_ptr< mfem::Coefficient > temp_bdr_coef_
 Coefficient containing the essential boundary values.
 
mfem_ext::StdFunctionOperator residual_with_bcs_
 mfem::Operator that describes the weight residual and its gradient with respect to temperature
 
std::unique_ptr< EquationSolvernonlin_solver_
 the specific methods and tolerances specified to solve the nonlinear residual equations
 
mfem_ext::FirstOrderODE ode_
 the ordinary differential equation that describes how to solve for the time derivative of temperature, given the current temperature and source terms
 
std::unique_ptr< mfem::HypreParMatrix > J_
 Assembled sparse matrix for the Jacobian.
 
std::unique_ptr< mfem::HypreParMatrix > J_e_
 
double dt_
 The current timestep.
 
double previous_dt_
 The previous timestep.
 
mfem::Vector u_
 Predicted temperature true dofs.
 
mfem::Vector u_predicted_
 Predicted temperature true dofs.
 
std::array< std::function< decltype((*residual_)(DifferentiateWRT< 1 >{}, 0.0, shape_displacement_, temperature_, temperature_rate_, *parameters_[parameter_indices].state...))(double)>, sizeof...(parameter_indices)> d_residual_d_
 Array functions computing the derivative of the residual with respect to each given parameter. More...
 
- Protected Attributes inherited from serac::BasePhysics
std::string name_ = {}
 Name of the physics module.
 
std::string mesh_tag_ = {}
 ID of the corresponding MFEMSidreDataCollection (denoting a mesh)
 
mfem::ParMesh & mesh_
 The primary mesh.
 
MPI_Comm comm_
 The MPI communicator.
 
std::vector< const serac::FiniteElementState * > states_
 List of finite element primal states associated with this physics module.
 
std::vector< const serac::FiniteElementState * > adjoints_
 List of finite element adjoint states associated with this physics module.
 
std::vector< const serac::FiniteElementDual * > duals_
 List of finite element duals associated with this physics module.
 
std::vector< ParameterInfoparameters_
 A vector of the parameters associated with this physics module.
 
FiniteElementStateshape_displacement_
 The parameter info associated with the shape displacement field. More...
 
std::unique_ptr< FiniteElementDualshape_displacement_sensitivity_
 Sensitivity with respect to the shape displacement field. More...
 
std::unordered_map< std::string, std::vector< serac::FiniteElementState > > checkpoint_states_
 A map containing optionally in-memory checkpointed primal states for transient adjoint solvers.
 
std::unordered_map< std::string, serac::FiniteElementStatecached_checkpoint_states_
 A container relating a checkpointed cycle and the associated finite element state fields. More...
 
std::optional< int > cached_checkpoint_cycle_
 An optional int for disk-based checkpointing containing the cycle number of the last retrieved checkpoint.
 
bool is_quasistatic_ = true
 Whether the simulation is time-independent.
 
double time_
 Current time for the forward pass.
 
double max_time_
 The maximum time reached for the forward solver.
 
double min_time_
 The time the forward solver was initialized to.
 
std::vector< double > timesteps_
 A vector of the timestep sizes (i.e. \(\Delta t\)) taken by the forward solver.
 
int cycle_
 Current cycle (forward pass time iteration count)
 
int max_cycle_
 The maximum cycle (forward pass iteration count) reached by the forward solver.
 
int min_cycle_
 The cycle the forward solver was initialized to.
 
double ode_time_point_
 The value of time at which the ODE solver wants to evaluate the residual.
 
int mpi_rank_
 MPI rank.
 
int mpi_size_
 MPI size.
 
std::unique_ptr< mfem::ParaViewDataCollection > paraview_dc_
 DataCollection pointer for optional paraview output.
 
std::unordered_map< std::string, std::unique_ptr< mfem::ParGridFunction > > paraview_dual_grid_functions_
 A optional map of the dual names and duals in grid function form for paraview output.
 
std::unique_ptr< mfem::ParGridFunction > shape_sensitivity_grid_function_
 A optional view of the shape sensitivity in grid function form for paraview output.
 
BoundaryConditionManager bcs_
 Boundary condition manager instance.
 
bool checkpoint_to_disk_
 A flag denoting whether to save the state to disk or memory as needed for dynamic adjoint solves.
 

Additional Inherited Members

- Protected Member Functions inherited from serac::BasePhysics
void CreateParaviewDataCollection () const
 Create a paraview data collection for the physics package if requested.
 
void UpdateParaviewDataCollection (const std::string &paraview_output_dir) const
 Update the paraview states, duals, parameters, and metadata (cycle, time) in preparation for output. More...
 
void initializeBasePhysicsStates (int cycle, double time)
 Protected, non-virtual method to reset physics states to zero. This does not reset design parameters or shape. More...
 
- Static Protected Attributes inherited from serac::BasePhysics
static constexpr int FLOAT_PRECISION_ = 8
 Number of significant figures to output for floating-point.
 

Detailed Description

template<int order, int dim, typename... parameter_space, int... parameter_indices>
class serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 87 of file heat_transfer.hpp.

Constructor & Destructor Documentation

◆ HeatTransfer() [1/3]

template<int order, int dim, typename... parameter_space, int... parameter_indices>
serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::HeatTransfer ( const NonlinearSolverOptions  nonlinear_opts,
const LinearSolverOptions  lin_opts,
const serac::TimesteppingOptions  timestepping_opts,
const std::string &  physics_name,
std::string  mesh_tag,
std::vector< std::string >  parameter_names = {},
int  cycle = 0,
double  time = 0.0,
bool  checkpoint_to_disk = false 
)
inline

Construct a new heat transfer object.

Parameters
[in]nonlinear_optsThe nonlinear solver options for solving the nonlinear residual equations
[in]lin_optsThe linear solver options for solving the linearized Jacobian equations
[in]timestepping_optsThe timestepping options for the heat transfer ordinary differential equations
[in]physics_nameA name for the physics module instance
[in]mesh_tagThe tag for the mesh in the StateManager to construct the physics module on
[in]parameter_namesA vector of the names of the requested parameter fields
[in]cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
[in]timeThe simulation time to initialize the physics module to
[in]checkpoint_to_diskA flag to save the transient states on disk instead of memory for the transient adjoint solves
Note
On parallel file systems (e.g. lustre), significant slowdowns and occasional errors were observed when writing and reading the needed trainsient states to disk for adjoint solves

Definition at line 117 of file heat_transfer.hpp.

◆ HeatTransfer() [2/3]

template<int order, int dim, typename... parameter_space, int... parameter_indices>
serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::HeatTransfer ( std::unique_ptr< serac::EquationSolver solver,
const serac::TimesteppingOptions  timestepping_opts,
const std::string &  physics_name,
std::string  mesh_tag,
std::vector< std::string >  parameter_names = {},
int  cycle = 0,
double  time = 0.0,
bool  checkpoint_to_disk = false 
)
inline

Construct a new heat transfer object.

Parameters
[in]solverThe nonlinear equation solver for the heat transfer equations
[in]timestepping_optsThe timestepping options for the heat transfer ordinary differential equations
[in]physics_nameA name for the physics module instance
[in]mesh_tagThe tag for the mesh in the StateManager to construct the physics module on
[in]parameter_namesA vector of the names of the requested parameter fields
[in]cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
[in]timeThe simulation time to initialize the physics module to
[in]checkpoint_to_diskA flag to save the transient states on disk instead of memory for the transient adjoint solves
Note
On parallel file systems (e.g. lustre), significant slowdowns and occasional errors were observed when writing and reading the needed trainsient states to disk for adjoint solves

Definition at line 142 of file heat_transfer.hpp.

◆ HeatTransfer() [3/3]

template<int order, int dim, typename... parameter_space, int... parameter_indices>
serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::HeatTransfer ( const HeatTransferInputOptions input_options,
const std::string &  physics_name,
const std::string &  mesh_tag,
int  cycle = 0,
double  time = 0.0 
)
inline

Construct a new Nonlinear HeatTransfer Solver object.

Parameters
[in]input_optionsThe solver information parsed from the input file
[in]physics_nameA name for the physics module instance
[in]mesh_tagThe tag for the mesh in the StateManager to construct the physics module on
[in]cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
[in]timeThe simulation time to initialize the physics module to

Definition at line 231 of file heat_transfer.hpp.

Member Function Documentation

◆ addCustomBoundaryIntegral()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
template<int... active_parameters, typename callable >
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::addCustomBoundaryIntegral ( DependsOn< active_parameters... >  ,
callable  qfunction,
const std::optional< Domain > &  optional_domain = std::nullopt 
)
inline

register a custom boundary integral calculation as part of the residual

Template Parameters
active_parametersa list of indices, describing which parameters to pass to the q-function
Parameters
qfunctiona callable that returns the normal heat flux on a boundary surface
optional_domainThe domain over which the integral is computed
heat_transfer.addCustomBoundaryIntegral(
DependsOn<>{},
[](double t, auto position, auto temperature, auto temperature_rate) {
auto [T, dT_dxi] = temperature;
auto q = 5.0*(T-25.0);
return q; // define a temperature-proportional heat-flux
});
Note
This method must be called prior to completeSetup()

Definition at line 680 of file heat_transfer.hpp.

◆ adjoint()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
const FiniteElementState& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::adjoint ( const std::string &  state_name) const
inlineoverridevirtual

Accessor for getting named finite element state adjoint solution from the physics modules.

Parameters
state_nameThe name of the Finite Element State adjoint solution to retrieve
Returns
The named adjoint Finite Element State

Implements serac::BasePhysics.

Definition at line 695 of file heat_transfer.hpp.

◆ advanceTimestep()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::advanceTimestep ( double  dt)
inlineoverridevirtual

Advance the heat conduction physics module in time.

Advance the underlying ODE with the requested time integration scheme using the previously set timestep.

Parameters
dtThe increment of simulation time to advance the underlying heat transfer problem

Implements serac::BasePhysics.

Definition at line 335 of file heat_transfer.hpp.

◆ completeSetup()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::completeSetup ( )
inlineoverridevirtual

Complete the initialization and allocation of the data structures.

This must be called before advanceTimestep().

Implements serac::BasePhysics.

Definition at line 711 of file heat_transfer.hpp.

◆ computeInitialConditionSensitivity()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
const std::unordered_map<std::string, const serac::FiniteElementDual&> serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::computeInitialConditionSensitivity ( )
inlineoverridevirtual

Compute the implicit sensitivity of the quantity of interest with respect to the initial temperature.

Returns
The sensitivity with respect to the initial temperature
Precondition
reverseAdjointTimestep must be called as many times as the forward solver was advanced before this is called

Reimplemented from serac::BasePhysics.

Definition at line 1014 of file heat_transfer.hpp.

◆ computeTimestepSensitivity()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
FiniteElementDual& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::computeTimestepSensitivity ( size_t  parameter_field)
inlineoverridevirtual

Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoint problem with respect to the parameter field for the last computed adjoint timestep.

Template Parameters
parameter_fieldThe index of the parameter to take a derivative with respect to
Returns
The sensitivity with respect to the parameter
Precondition
reverseAdjointTimestep with an appropriate adjoint load must be called prior to this method.

Reimplemented from serac::BasePhysics.

Definition at line 974 of file heat_transfer.hpp.

◆ computeTimestepShapeSensitivity()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
FiniteElementDual& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::computeTimestepShapeSensitivity ( )
inlineoverridevirtual

Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoint problem with respect to the shape displacement field for the last computed adjoint timestep.

Returns
The sensitivity with respect to the shape displacement
Precondition
reverseAdjointTimestep with an appropriate adjoint load must be called prior to this method.

Reimplemented from serac::BasePhysics.

Definition at line 994 of file heat_transfer.hpp.

◆ getCheckpointedStates()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
std::unordered_map<std::string, FiniteElementState> serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::getCheckpointedStates ( int  cycle_to_load) const
inlineoverridevirtual

Accessor for getting named finite element state primal solution from the physics modules at a given checkpointed cycle index.

Parameters
cycle_to_loadThe previous timestep where the state solution is requested
Returns
The named primal Finite Element State

Reimplemented from serac::BasePhysics.

Definition at line 946 of file heat_transfer.hpp.

◆ initializeThermalStates()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::initializeThermalStates ( )
inline

Non virtual method to reset thermal states to zero. This does not reset design parameters or shape.

Parameters
[in]cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
[in]timeThe simulation time to initialize the physics module to

Definition at line 279 of file heat_transfer.hpp.

◆ resetStates()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::resetStates ( int  cycle = 0,
double  time = 0.0 
)
inlineoverridevirtual

Method to reset physics states to zero. This does not reset design parameters or shape.

Parameters
[in]cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
[in]timeThe simulation time to initialize the physics module to

Implements serac::BasePhysics.

Definition at line 306 of file heat_transfer.hpp.

◆ reverseAdjointTimestep()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::reverseAdjointTimestep ( )
inlineoverridevirtual

Solve the adjoint problem.

Precondition
It is expected that the forward analysis is complete and the current temperature state is valid
It is expected that the adjoint load has already been set in HeatTransfer::setAdjointLoad
Note
If the essential boundary dual is not specified, homogeneous essential boundary conditions are applied to the adjoint system
We provide a quick derivation for the discrete adjoint equations and notation used here for backward Euler. There are two equations satisfied at each step of the forward solve using backward Euler
1). \(r^n(u^n, v^n; d) = 0,\)
where \(u^n\) is the end-step primal value, \(d\) are design parameters, and \(v^n\) is the central difference velocity, satisfying
2). \(\Delta t \, v^n = u^n-u^{n-1}.\)
We are interesting in the implicit sensitivity of a qoi (quantity of interest), \(\sum_{n=1}^N \pi^n(u^n, v^n; d)\), while maintaining the above constraints. We construct a Lagrangian that adds 'zero' to this qoi using the free multipliers \(\lambda^n\) (which we call the adjoint temperature) and \(\mu^n\) (which can eventually be intepreted as the implicit sensitivity of the qoi with respect to the start-of-step primal value \(u^{n-1}\)) with

\[ \mathcal{L} := \sum_{n=1}^N \pi(u^n, v^n; d) + \lambda^n \cdot r^n(u^n, v^n; d) + \mu^n \cdot (\Delta t \, v^n - u^n + u^{n-1}).\]

We are interesting in the total derivative

\[\frac{d\mathcal{L}}{dp} = \sum_{n=1}^N \pi^n_{,u} u^n_{,d} + \pi^n_{,v} v^n_{,d} + \pi^n_{,d} + \lambda^n \cdot \left( r^n_{,u} u^n_{,d} + r^n_{,v} v^n_{,d} + r^n_{,d} \right) + \mu^n \cdot (\Delta t \, v^n_{,d} - u^n_{,d} + u^{n-1}_{,d}). \]

We are free to choose \(\lambda^n\) and \(\mu^n\) in any convenient way, and the way we choose is by grouping terms involving \(u^n_{,d}\) and \(v^n_{,d}\), and setting the multipliers such that those terms cancel out in the final expression of the qoi sensitivity. In particular, by choosing

\[ \lambda^n = - \left[ r_{,u}^n + \frac{r_{,v}^n}{\Delta t} \right]^{-T} \left( \pi^n_{,u} + \frac{\pi^n_{,v}}{\Delta t} + \mu^{n+1} \right) \]

and

\[ \mu^n = -\frac{1}{\Delta t}\left( \pi^n_{,v} + \lambda^n \cdot r^n_{,v} \right), \]

we find

\[ \frac{d\mathcal{L}}{dp} = \sum_{n=1}^N \left( \pi^n_{,d} + \lambda^n \cdot r^n_{,d} \right) + \mu^1 \cdot u^{0}_{,d}, \]

where the multiplier/adjoint equations are solved backward in time starting at \(n=N\), \(\mu^{N+1} = 0\), and \(u^{0}_{,d}\) is the sensitivity of the initial primal variable with respect to the parameters.
We call the quantities \(\pi^n_{,u}\) and \(\pi^n_{,v}\) the temperature and temperature-rate adjoint loads, respectively.

Reimplemented from serac::BasePhysics.

Definition at line 851 of file heat_transfer.hpp.

◆ setAdjointLoad()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
virtual void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setAdjointLoad ( std::unordered_map< std::string, const serac::FiniteElementDual & >  loads)
inlineoverridevirtual

Set the loads for the adjoint reverse timestep solve.

Parameters
loadsThe loads (e.g. right hand sides) for the adjoint problem
Precondition
The adjoint load map is expected to contain an entry named "temperature"
The adjoint load map may contain an entry named "temperature_rate"

These loads are typically defined as derivatives of a downstream quantity of intrest with respect to a primal solution field (in this case, temperature). For this physics module, the unordered map is expected to have two entries with the keys "temperature" and "temperature_rate". Note that the "temperature_rate" load is only used by transient (i.e. non-quasi-static) adjoint calculations.

Reimplemented from serac::BasePhysics.

Definition at line 801 of file heat_transfer.hpp.

◆ setFluxBCs()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
template<int... active_parameters, typename FluxType >
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setFluxBCs ( DependsOn< active_parameters... >  ,
FluxType  flux_function,
const std::optional< Domain > &  optional_domain = std::nullopt 
)
inline

Set the thermal flux boundary condition.

Template Parameters
FluxTypeThe type of the thermal flux object
Parameters
flux_functionA function describing the flux applied to a boundary
optional_domainThe domain over which the flux is applied. If nothing is supplied the entire boundary is used.
Precondition
FluxType must be a object that can be called with the following arguments:
  1. tensor<T,dim> x the spatial coordinates for the quadrature point
  2. tensor<T,dim> n the outward-facing unit normal for the quadrature point
  3. double t the time (note: time will be handled differently in the future)
  4. T temperature the current temperature at the quadrature point
  1. tuple{value, derivative}, a variadic list of tuples (each with a values and derivative), one tuple for each of the trial spaces specified in the DependsOn<...> argument.
Note
The actual types of these arguments passed will be double, tensor<double, ... > or tuples thereof when doing direct evaluation. When differentiating with respect to one of the inputs, its stored values will change to dual numbers rather than double. (e.g. tensor<double,3> becomes tensor<dual<...>, 3>)
This method must be called prior to completeSetup()

Definition at line 547 of file heat_transfer.hpp.

◆ setMaterial()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
template<int... active_parameters, typename MaterialType >
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setMaterial ( DependsOn< active_parameters... >  ,
MaterialType  material 
)
inline

Set the thermal material model for the physics solver.

Template Parameters
MaterialTypeThe thermal material type
Parameters
materialA material containing heat capacity and thermal flux evaluation information
Precondition
material must be a object that can be called with the following arguments:
  1. tensor<T,dim> x the spatial position of the material evaluation call
  2. T temperature the current temperature at the quadrature point
  3. tensor<T,dim> the spatial gradient of the temperature at the quadrature point
  4. tuple{value, derivative}, a tuple of values and derivatives for each parameter field specified in the DependsOn<...> argument.
Note
The actual types of these arguments passed will be double, tensor<double, ... > or tuples thereof when doing direct evaluation. When differentiating with respect to one of the inputs, its stored values will change to dual numbers rather than double. (e.g. tensor<double,3> becomes tensor<dual<...>, 3>)
Precondition
MaterialType must return a serac::tuple of volumetric heat capacity and thermal flux when operator() is called with the arguments listed above.
Note
This method must be called prior to completeSetup()

Definition at line 439 of file heat_transfer.hpp.

◆ setShapeDisplacement()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setShapeDisplacement ( std::function< void(const mfem::Vector &x, mfem::Vector &shape_disp)>  shape_disp)
inline

Set the underlying finite element state to a prescribed shape displacement.

This field will perturb the mesh nodes as required by shape optimization workflows.

Parameters
shape_dispThe function describing the shape displacement field

Definition at line 577 of file heat_transfer.hpp.

◆ setSource()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
template<int... active_parameters, typename SourceType >
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setSource ( DependsOn< active_parameters... >  ,
SourceType  source_function,
const std::optional< Domain > &  optional_domain = std::nullopt 
)
inline

Set the thermal source function.

Template Parameters
SourceTypeThe type of the source function
Parameters
source_functionA source function for a prescribed thermal load
optional_domainThe domain over which the source is applied. If nothing is supplied the entire domain is used.
Precondition
source_function must be a object that can be called with the following arguments:
  1. tensor<T,dim> x the spatial coordinates for the quadrature point
  2. double t the time (note: time will be handled differently in the future)
  3. T temperature the current temperature at the quadrature point
  4. tensor<T,dim> the spatial gradient of the temperature at the quadrature point
  5. tuple{value, derivative}, a variadic list of tuples (each with a values and derivative), one tuple for each of the trial spaces specified in the DependsOn<...> argument.
Note
The actual types of these arguments passed will be double, tensor<double, ... > or tuples thereof when doing direct evaluation. When differentiating with respect to one of the inputs, its stored values will change to dual numbers rather than double. (e.g. tensor<double,3> becomes tensor<dual<...>, 3>)
This method must be called prior to completeSetup()

Definition at line 497 of file heat_transfer.hpp.

◆ setState()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setState ( const std::string &  state_name,
const FiniteElementState state 
)
inlineoverridevirtual

Set the primal solution field (temperature) for the underlying heat transfer solver.

Parameters
state_nameThe name of the field to initialize (must be "temperature")
stateThe finite element state vector containing the values for either the temperature field

It is expected that state has the same underlying finite element space and mesh as the selected primal solution field.

Implements serac::BasePhysics.

Definition at line 629 of file heat_transfer.hpp.

◆ setTemperature()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setTemperature ( std::function< double(const mfem::Vector &x, double t)>  temp)
inline

Set the underlying finite element state to a prescribed temperature.

Parameters
tempThe function describing the temperature field
Note
This will override any existing solution values in the temperature field

Definition at line 461 of file heat_transfer.hpp.

◆ setTemperatureBCs()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
void serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::setTemperatureBCs ( const std::set< int > &  temp_bdr,
std::function< double(const mfem::Vector &x, double t)>  temp 
)
inline

Set essential temperature boundary conditions (strongly enforced)

Parameters
[in]temp_bdrThe boundary attributes on which to enforce a temperature
[in]tempThe prescribed boundary temperature function
Note
This should be called prior to completeSetup()

Definition at line 320 of file heat_transfer.hpp.

◆ state()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
const FiniteElementState& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::state ( const std::string &  state_name) const
inlineoverridevirtual

Accessor for getting named finite element state primal solution from the physics modules.

Parameters
state_nameThe name of the Finite Element State primal solution to retrieve
Returns
The named primal Finite Element State

Implements serac::BasePhysics.

Definition at line 607 of file heat_transfer.hpp.

◆ stateNames()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
virtual std::vector<std::string> serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::stateNames ( ) const
inlineoverridevirtual

Get a vector of the finite element state primal solution names.

Returns
The primal solution names

Implements serac::BasePhysics.

Definition at line 649 of file heat_transfer.hpp.

◆ temperature()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
const serac::FiniteElementState& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::temperature ( ) const
inline

Get the temperature state.

Returns
A reference to the current temperature finite element state

Definition at line 592 of file heat_transfer.hpp.

◆ temperatureRate()

template<int order, int dim, typename... parameter_space, int... parameter_indices>
const serac::FiniteElementState& serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::temperatureRate ( ) const
inline

Get the temperature rate of change state.

Returns
A reference to the current temperature rate of change finite element state

Definition at line 599 of file heat_transfer.hpp.

Member Data Documentation

◆ d_residual_d_

template<int order, int dim, typename... parameter_space, int... parameter_indices>
std::array<std::function<decltype((*residual_)(DifferentiateWRT<1>{}, 0.0, shape_displacement_, temperature_, temperature_rate_, *parameters_[parameter_indices].state...))(double)>, sizeof...(parameter_indices)> serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::d_residual_d_
protected
Initial value:
= {[&](double TIME) {
return (*residual_)(DifferentiateWRT<NUM_STATE_VARS + 1 + parameter_indices>{}, TIME, shape_displacement_,
temperature_, temperature_rate_, *parameters_[parameter_indices].state...);
}...}
FiniteElementState & shape_displacement_
The parameter info associated with the shape displacement field.
std::vector< ParameterInfo > parameters_
A vector of the parameters associated with this physics module.
std::unique_ptr< ShapeAwareFunctional< shape_trial, test(scalar_trial, scalar_trial, parameter_space...)> > residual_
serac::Functional that is used to calculate the residual and its derivatives
FiniteElementState temperature_rate_
Rate of change in temperature at the current adjoint timestep.

Array functions computing the derivative of the residual with respect to each given parameter.

Note
This is needed so the user can ask for a specific sensitivity at runtime as opposed to it being a template parameter.

Definition at line 1100 of file heat_transfer.hpp.

◆ J_e_

template<int order, int dim, typename... parameter_space, int... parameter_indices>
std::unique_ptr<mfem::HypreParMatrix> serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >::J_e_
protected

rows and columns of J_ that have been separated out because are associated with essential boundary conditions

Definition at line 1080 of file heat_transfer.hpp.


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