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

The operator-split thermal-structural solver. More...

#include <thermomechanics.hpp>

Collaboration diagram for serac::Thermomechanics< order, dim, parameter_space >:
Collaboration graph
[legend]

Classes

struct  MechanicalMaterialInterface
 This is an adaptor class that makes a thermomechanical material usable by the solid mechanics module, by discarding the thermal-specific information. More...
 
struct  ThermalMaterialInterface
 This is an adaptor class that makes a thermomechanical material usable by the thermal module, by discarding the solid-mechanics-specific information. More...
 

Public Member Functions

 Thermomechanics (const NonlinearSolverOptions thermal_nonlin_opts, const LinearSolverOptions thermal_lin_opts, TimesteppingOptions thermal_timestepping, const NonlinearSolverOptions solid_nonlin_opts, const LinearSolverOptions solid_lin_opts, TimesteppingOptions solid_timestepping, GeometricNonlinearities geom_nonlin, const std::string &physics_name, std::string mesh_tag, int cycle=0, double time=0.0)
 Construct a new coupled Thermal-SolidMechanics object. More...
 
 Thermomechanics (std::unique_ptr< EquationSolver > thermal_solver, TimesteppingOptions thermal_timestepping, std::unique_ptr< EquationSolver > solid_solver, TimesteppingOptions solid_timestepping, GeometricNonlinearities geom_nonlin, const std::string &physics_name, std::string mesh_tag, int cycle=0, double time=0.0)
 Construct a new coupled Thermal-SolidMechanics object. More...
 
 Thermomechanics (const HeatTransferInputOptions &thermal_options, const SolidMechanicsInputOptions &solid_options, const std::string &physics_name, std::string mesh_tag, int cycle=0, double time=0.0)
 Construct a new Thermal-SolidMechanics Functional object from input file options. More...
 
 Thermomechanics (const ThermomechanicsInputOptions &options, const std::string &physics_name, std::string mesh_tag, int cycle=0, double time=0.0)
 Construct a new Thermal-SolidMechanics Functional object from input file options. More...
 
void completeSetup () override
 Complete the initialization and allocation of the data structures. 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...
 
const FiniteElementStatestate (const std::string &state_name) const override
 Accessor for getting named finite element state primal fields from the physics modules. More...
 
void setState (const std::string &state_name, const FiniteElementState &state) override
 Set the primal solution field (displacement, velocity, temperature) for the underlying thermomechanics solver. More...
 
virtual std::vector< std::string > stateNames () const override
 Get a vector of the finite element state solution variable names. More...
 
const FiniteElementStateadjoint (const std::string &state_name) const override
 Accessor for getting named finite element adjoint fields from the physics modules. More...
 
void advanceTimestep (double dt) override
 Advance the timestep. More...
 
template<typename T >
std::shared_ptr< QuadratureData< T > > createQuadratureDataBuffer (T initial_state)
 Create a shared ptr to a quadrature data buffer for the given material type. More...
 
template<int... active_parameters, typename MaterialType , typename StateType >
void setMaterial (DependsOn< active_parameters... >, MaterialType material, std::shared_ptr< QuadratureData< StateType >> qdata)
 Set the thermomechanical material response. More...
 
template<typename MaterialType , typename StateType = Empty>
void setMaterial (MaterialType material, std::shared_ptr< QuadratureData< StateType >> qdata=EmptyQData)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
void setTemperatureBCs (const std::set< int > &temperature_attributes, std::function< double(const mfem::Vector &x, double t)> prescribed_value)
 Set essential temperature boundary conditions (strongly enforced) More...
 
void setDisplacementBCs (const std::set< int > &displacement_attributes, std::function< void(const mfem::Vector &x, mfem::Vector &disp)> prescribed_value)
 Set essential displacement boundary conditions (strongly enforced) More...
 
template<typename FluxType >
void setHeatFluxBCs (FluxType flux_function)
 Set the thermal flux boundary condition. More...
 
void setDisplacement (std::function< void(const mfem::Vector &x, mfem::Vector &u)> displacement)
 Set the underlying finite element state to a prescribed displacement. More...
 
void setTemperature (std::function< double(const mfem::Vector &x, double t)> temperature)
 Set the underlying finite element state to a prescribed temperature. More...
 
template<typename BodyForceType >
void addBodyForce (BodyForceType body_force_function)
 Set the body forcefunction. More...
 
template<typename HeatSourceType >
void addHeatSource (HeatSourceType source_function)
 Set the thermal source function. More...
 
const serac::FiniteElementStatedisplacement () const
 Get the displacement state. More...
 
const serac::FiniteElementStatetemperature () const
 Get the temperature state. More...
 
- 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 const FiniteElementDualcomputeTimestepSensitivity (size_t)
 Compute the implicit sensitivity of the quantity of interest used in defining the adjoint load with respect to the parameter field (d QOI/d state * d state/d parameter). More...
 
virtual const FiniteElementDualcomputeTimestepShapeSensitivity ()
 Compute the implicit sensitivity of the quantity of interest used in defining the adjoint load with respect to the shape displacement field (d QOI/d state * d state/d shape displacement). More...
 
virtual const std::unordered_map< std::string, const serac::FiniteElementDual & > computeInitialConditionSensitivity ()
 Compute the implicit sensitivity of the quantity of interest with respect to the initial condition fields. More...
 
virtual void setAdjointLoad (std::unordered_map< std::string, const serac::FiniteElementDual & >)
 Set the loads for the adjoint reverse timestep solve.
 
virtual void reverseAdjointTimestep ()
 Solve the adjoint reverse timestep problem. 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.
 

Protected Types

using displacement_field = H1< order, dim >
 the function space for the displacement field
 
using temperature_field = H1< order >
 the function space for the temperature field
 

Protected Attributes

HeatTransfer< order, dim, Parameters< displacement_field, parameter_space... > > thermal_
 Submodule to compute the heat transfer physics.
 
SolidMechanics< order, dim, Parameters< temperature_field, parameter_space... > > solid_
 Submodule to compute the mechanics.
 
- 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...
 
virtual std::unordered_map< std::string, FiniteElementStategetCheckpointedStates (int cycle) const
 Accessor for getting all of the primal solutions from the physics modules at a given checkpointed cycle index. 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>
class serac::Thermomechanics< order, dim, parameter_space >

The operator-split thermal-structural solver.

Uses Functional to compute action of operators

Definition at line 32 of file thermomechanics.hpp.

Constructor & Destructor Documentation

◆ Thermomechanics() [1/4]

template<int order, int dim, typename... parameter_space>
serac::Thermomechanics< order, dim, parameter_space >::Thermomechanics ( const NonlinearSolverOptions  thermal_nonlin_opts,
const LinearSolverOptions  thermal_lin_opts,
TimesteppingOptions  thermal_timestepping,
const NonlinearSolverOptions  solid_nonlin_opts,
const LinearSolverOptions  solid_lin_opts,
TimesteppingOptions  solid_timestepping,
GeometricNonlinearities  geom_nonlin,
const std::string &  physics_name,
std::string  mesh_tag,
int  cycle = 0,
double  time = 0.0 
)
inline

Construct a new coupled Thermal-SolidMechanics object.

Parameters
thermal_nonlin_optsThe options for solving the nonlinear heat conduction residual equations
thermal_lin_optsThe options for solving the linearized Jacobian heat transfer equations
thermal_timesteppingThe timestepping options for the heat transfer operator
solid_nonlin_optsThe options for solving the nonlinear solid mechanics residual equations
solid_lin_optsThe options for solving the linearized Jacobian solid mechanics equations
solid_timesteppingThe timestepping options for the solid solver
geom_nonlinFlag to include geometric nonlinearities
physics_nameA name for the physics module instance
mesh_tagThe tag for the mesh in the StateManager to construct the physics module on
cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
timeThe simulation time to initialize the physics module to

Definition at line 49 of file thermomechanics.hpp.

◆ Thermomechanics() [2/4]

template<int order, int dim, typename... parameter_space>
serac::Thermomechanics< order, dim, parameter_space >::Thermomechanics ( std::unique_ptr< EquationSolver thermal_solver,
TimesteppingOptions  thermal_timestepping,
std::unique_ptr< EquationSolver solid_solver,
TimesteppingOptions  solid_timestepping,
GeometricNonlinearities  geom_nonlin,
const std::string &  physics_name,
std::string  mesh_tag,
int  cycle = 0,
double  time = 0.0 
)
inline

Construct a new coupled Thermal-SolidMechanics object.

Parameters
thermal_solverThe nonlinear equation solver for the heat conduction equations
thermal_timesteppingThe timestepping options for the thermal solver
solid_solverThe nonlinear equation solver for the solid mechanics equations
solid_timesteppingThe timestepping options for the solid solver
geom_nonlinFlag to include geometric nonlinearities
physics_nameA name for the physics module instance
mesh_tagThe tag for the mesh in the StateManager to construct the physics module on
cycleThe simulation cycle (i.e. timestep iteration) to intialize the physics module to
timeThe simulation time to initialize the physics module to

Definition at line 76 of file thermomechanics.hpp.

◆ Thermomechanics() [3/4]

template<int order, int dim, typename... parameter_space>
serac::Thermomechanics< order, dim, parameter_space >::Thermomechanics ( const HeatTransferInputOptions thermal_options,
const SolidMechanicsInputOptions solid_options,
const std::string &  physics_name,
std::string  mesh_tag,
int  cycle = 0,
double  time = 0.0 
)
inline

Construct a new Thermal-SolidMechanics Functional object from input file options.

Parameters
[in]thermal_optionsThe thermal physics module input file option struct
[in]solid_optionsThe solid physics module input file option struct
[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 104 of file thermomechanics.hpp.

◆ Thermomechanics() [4/4]

template<int order, int dim, typename... parameter_space>
serac::Thermomechanics< order, dim, parameter_space >::Thermomechanics ( const ThermomechanicsInputOptions options,
const std::string &  physics_name,
std::string  mesh_tag,
int  cycle = 0,
double  time = 0.0 
)
inline

Construct a new Thermal-SolidMechanics Functional object from input file options.

Parameters
[in]optionsThe thermal solid physics module input file option struct
[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 122 of file thermomechanics.hpp.

Member Function Documentation

◆ addBodyForce()

template<int order, int dim, typename... parameter_space>
template<typename BodyForceType >
void serac::Thermomechanics< order, dim, parameter_space >::addBodyForce ( BodyForceType  body_force_function)
inline

Set the body forcefunction.

Template Parameters
BodyForceTypeThe type of the body force load
Precondition
body_force_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. 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>)

Definition at line 485 of file thermomechanics.hpp.

◆ addHeatSource()

template<int order, int dim, typename... parameter_space>
template<typename HeatSourceType >
void serac::Thermomechanics< order, dim, parameter_space >::addHeatSource ( HeatSourceType  source_function)
inline

Set the thermal source function.

Template Parameters
HeatSourceTypeThe type of the source function
Parameters
source_functionA source function for a prescribed thermal load
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>)

Definition at line 510 of file thermomechanics.hpp.

◆ adjoint()

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

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

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

Implements serac::BasePhysics.

Definition at line 223 of file thermomechanics.hpp.

◆ advanceTimestep()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::advanceTimestep ( double  dt)
inlineoverridevirtual

Advance the timestep.

Parameters
dtThe increment of simulation time to advance the underlying thermomechanical problem

Implements serac::BasePhysics.

Definition at line 241 of file thermomechanics.hpp.

◆ completeSetup()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::completeSetup ( )
inlineoverridevirtual

Complete the initialization and allocation of the data structures.

Note
This must be called before AdvanceTimestep().

Implements serac::BasePhysics.

Definition at line 139 of file thermomechanics.hpp.

◆ createQuadratureDataBuffer()

template<int order, int dim, typename... parameter_space>
template<typename T >
std::shared_ptr<QuadratureData<T> > serac::Thermomechanics< order, dim, parameter_space >::createQuadratureDataBuffer ( initial_state)
inline

Create a shared ptr to a quadrature data buffer for the given material type.

Template Parameters
Tthe type to be created at each quadrature point
Parameters
initial_statethe value to be broadcast to each quadrature point
Returns
std::shared_ptr< QuadratureData<T> >

Definition at line 261 of file thermomechanics.hpp.

◆ displacement()

template<int order, int dim, typename... parameter_space>
const serac::FiniteElementState& serac::Thermomechanics< order, dim, parameter_space >::displacement ( ) const
inline

Get the displacement state.

Returns
A reference to the current displacement finite element state

Definition at line 520 of file thermomechanics.hpp.

◆ resetStates()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::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 151 of file thermomechanics.hpp.

◆ setDisplacement()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::setDisplacement ( std::function< void(const mfem::Vector &x, mfem::Vector &u)>  displacement)
inline

Set the underlying finite element state to a prescribed displacement.

Parameters
displacementThe function describing the displacement field

Definition at line 454 of file thermomechanics.hpp.

◆ setDisplacementBCs()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::setDisplacementBCs ( const std::set< int > &  displacement_attributes,
std::function< void(const mfem::Vector &x, mfem::Vector &disp)>  prescribed_value 
)
inline

Set essential displacement boundary conditions (strongly enforced)

Parameters
[in]displacement_attributesThe boundary attributes on which to enforce a displacement
[in]prescribed_valueThe prescribed boundary displacement function

Definition at line 414 of file thermomechanics.hpp.

◆ setHeatFluxBCs()

template<int order, int dim, typename... parameter_space>
template<typename FluxType >
void serac::Thermomechanics< order, dim, parameter_space >::setHeatFluxBCs ( FluxType  flux_function)
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
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>)
: until mfem::GetFaceGeometricFactors implements their JACOBIANS option, (or we implement a replacement kernel ourselves) we are not able to compute shape sensitivities for boundary integrals.

Definition at line 444 of file thermomechanics.hpp.

◆ setMaterial()

template<int order, int dim, typename... parameter_space>
template<int... active_parameters, typename MaterialType , typename StateType >
void serac::Thermomechanics< order, dim, parameter_space >::setMaterial ( DependsOn< active_parameters... >  ,
MaterialType  material,
std::shared_ptr< QuadratureData< StateType >>  qdata 
)
inline

Set the thermomechanical material response.

Template Parameters
MaterialTypeThe thermomechanical material type
StateTypeThe type that contains the internal variables for MaterialType
Parameters
materialA material that provides a function to evaluate stress, heat flux, density, and heat capacity
qdatathe buffer of material internal variables at each quadrature point
Precondition
material must be a object that can be called with the following arguments:
  1. MaterialType::State & state an mutable reference to the internal variables for this quadrature point
  2. tensor<T,dim,dim> du_dx the displacement gradient at this quadrature point
  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 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 Cauchy stress, volumetric heat capacity, internal heat source, and thermal flux when operator() is called with the arguments listed above.

Definition at line 378 of file thermomechanics.hpp.

◆ setState()

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

Set the primal solution field (displacement, velocity, temperature) for the underlying thermomechanics solver.

Parameters
state_nameThe name of the field to initialize ("displacement", or "velocity")
stateThe finite element state vector containing the values for either the displacement or velocity fields

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 189 of file thermomechanics.hpp.

◆ setTemperature()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::setTemperature ( std::function< double(const mfem::Vector &x, double t)>  temperature)
inline

Set the underlying finite element state to a prescribed temperature.

Parameters
temperatureThe function describing the temperature field

Definition at line 464 of file thermomechanics.hpp.

◆ setTemperatureBCs()

template<int order, int dim, typename... parameter_space>
void serac::Thermomechanics< order, dim, parameter_space >::setTemperatureBCs ( const std::set< int > &  temperature_attributes,
std::function< double(const mfem::Vector &x, double t)>  prescribed_value 
)
inline

Set essential temperature boundary conditions (strongly enforced)

Parameters
[in]temperature_attributesThe boundary attributes on which to enforce a temperature
[in]prescribed_valueThe prescribed boundary temperature function

Definition at line 402 of file thermomechanics.hpp.

◆ state()

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

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

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

Implements serac::BasePhysics.

Definition at line 164 of file thermomechanics.hpp.

◆ stateNames()

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

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

Returns
The solution variable names

Implements serac::BasePhysics.

Definition at line 212 of file thermomechanics.hpp.

◆ temperature()

template<int order, int dim, typename... parameter_space>
const serac::FiniteElementState& serac::Thermomechanics< order, dim, parameter_space >::temperature ( ) const
inline

Get the temperature state.

Returns
A reference to the current temperature finite element state

Definition at line 527 of file thermomechanics.hpp.


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