|
Serac
0.1
Serac is an implicit thermal strucural mechanics simulation code.
|
This is the abstract base class for a generic forward solver. More...
#include <base_physics.hpp>

Classes | |
| struct | ParameterInfo |
| The information needed for the physics parameters stored as Finite Element State fields. More... | |
Public Member Functions | |
| BasePhysics (std::string physics_name, std::shared_ptr< serac::Mesh > mesh, 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 const std::vector< double > & | timesteps () const |
| Get a vector of the timestep sizes (i.e. \(\Delta t\)s) taken by the forward solver. More... | |
| virtual void | resetStates (int cycle=0, double time=0.0)=0 |
| Base method to reset physics states to the initial time. This does not reset design parameters or shape. More... | |
| virtual void | resetAdjointStates () |
| Base method to reset physics states back to the end of time to start adjoint calculations again. This does not reset design parameters or shape. | |
| virtual void | completeSetup ()=0 |
| Complete the setup and allocate the necessary data structures. More... | |
| virtual const FiniteElementState & | state (const std::string &state_name) const =0 |
| Accessor for getting named finite element state primal solution from the physics modules. More... | |
| virtual void | setState (const std::string &, const FiniteElementState &)=0 |
| Set the primal solution field values of the underlying physics solver. | |
| virtual std::vector< std::string > | stateNames () const =0 |
| Get a vector of the finite element state primal solution names. More... | |
| virtual const FiniteElementState & | adjoint (const std::string &adjoint_name) const =0 |
| Accessor for getting named finite element state adjoint solution from the physics modules. More... | |
| virtual std::vector< std::string > | adjointNames () const |
| Get a vector of the finite element state adjoint solution names. More... | |
| virtual std::vector< std::string > | dualNames () const |
| Get a vector of the finite element state dual (reaction) solution names. More... | |
| virtual const FiniteElementDual & | dual (const std::string &dual_name) const |
| Accessor for getting named finite element state dual (reaction) solution from the physics modules. More... | |
| virtual const FiniteElementState & | dualAdjoint (const std::string &dual_name) const |
| Accessor for getting named finite element state dual adjoint (reaction adjoint load) from the physics modules. More... | |
| virtual const FiniteElementState & | shapeDisplacement () const |
| Accessor for getting the shape displacement field from the physics modules. More... | |
| virtual const FiniteElementState & | parameter (const std::string ¶meter_name) const |
| Accessor for getting named finite element state parameter fields from the physics modules. More... | |
| virtual const FiniteElementState & | parameter (std::size_t parameter_index) const |
| Accessor for getting indexed finite element state parameter fields from the physics modules. More... | |
| virtual std::vector< std::string > | parameterNames () const |
| Get a vector of the finite element state parameter names. More... | |
| virtual void | setParameter (const size_t parameter_index, const FiniteElementState ¶meter_state) |
| Deep copy a parameter field into the internally-owned parameter used for simulations. More... | |
| virtual void | setShapeDisplacement (const FiniteElementState &shape_displacement) |
| Set the current shape displacement for the underlying mesh. More... | |
| virtual FiniteElementDual | computeTimestepSensitivity (size_t parameter_index) |
| 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 FiniteElementDual & | computeTimestepShapeSensitivity () |
| 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 () const |
| Compute the implicit sensitivity of the quantity of interest with respect to the initial condition fields. More... | |
| virtual void | advanceTimestep (double dt)=0 |
| Advance the state variables according to the chosen time integrator. More... | |
| virtual void | setAdjointLoad (std::unordered_map< std::string, const serac::FiniteElementDual & > string_to_dual) |
| Set the loads for the adjoint reverse timestep solve. More... | |
| virtual void | setDualAdjointBcs (std::unordered_map< std::string, const serac::FiniteElementState & > string_to_bc) |
| Set the dual loads (dirichlet values) for the adjoint reverse timestep solve This must be called after setAdjointLoad. More... | |
| virtual void | reverseAdjointTimestep () |
| Solve the adjoint reverse timestep problem. More... | |
| virtual void | initializationStep () |
| Initialize any fields nessary for before the first step of the time integration. | |
| virtual void | reverseAdjointInitializationStep () |
| Compute adjoint sensitivities back through initializationStep. | |
| 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... | |
| void | loadCheckpointedStatesFromDisk (int cycle) |
| load checkpointed states from disk into states array More... | |
| virtual FiniteElementState | loadCheckpointedState (const std::string &state_name, int cycle) |
| Accessor for getting a single named finite element state primal solution from the physics modules at a given checkpointed cycle index. More... | |
| virtual FiniteElementDual | loadCheckpointedDual ([[maybe_unused]] const std::string &state_name, [[maybe_unused]] int cycle) |
| Accessor for getting a single named finite element dual 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 serac::Mesh & | mesh () const |
| Returns a reference to the mesh object. | |
| const mfem::ParMesh & | mfemParMesh () const |
| Returns a reference to the mfem ParMesh object. | |
| mfem::ParMesh & | mfemParMesh () |
| Returns a reference to the mfem ParMesh object. | |
| std::string | name () const |
| Return the name of the physics. | |
Protected Member Functions | |
| const FiniteElementDual & | shapeDisplacementSensitivity () const |
| Internally used accessor for getting the shape displacement sensitivity from the physics modules. More... | |
| void | CreateParaviewDataCollection () const |
| Create a paraview data collection for the physics package if requested. | |
| void | UpdateParaviewDataCollection (const std::string ¶view_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... | |
| std::unordered_map< std::string, FiniteElementState > | getCheckpointedStates (int cycle) |
| Accessor for getting all of the primal solutions from the physics modules at a given checkpointed cycle index. More... | |
Protected Attributes | |
| std::string | name_ = {} |
| Name of the physics module. | |
| std::shared_ptr< serac::Mesh > | 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< const serac::FiniteElementState * > | dual_adjoints_ |
| List of adjoint finite element duals associated with this physics module. | |
| std::vector< ParameterInfo > | parameters_ |
| A vector of the parameters associated with this physics module. | |
| 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::FiniteElementState > | cached_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 | dt_ |
| Current time step. | |
| 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. | |
| serac::FiniteElementState | shape_displacement_ |
| The shape displacement field. | |
| serac::FiniteElementDual | shape_displacement_dual_ |
| The shape displacement field sensitivity. | |
| 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. | |
Static Protected Attributes | |
| static constexpr int | FLOAT_PRECISION_ = 8 |
| Number of significant figures to output for floating-point. | |
This is the abstract base class for a generic forward solver.
Definition at line 61 of file base_physics.hpp.
| serac::BasePhysics::BasePhysics | ( | std::string | physics_name, |
| std::shared_ptr< serac::Mesh > | mesh, | ||
| int | cycle = 0, |
||
| double | time = 0.0, |
||
| bool | checkpoint_to_disk = false |
||
| ) |
Empty constructor.
| [in] | physics_name | Name of the physics module instance |
| [in] | mesh | The mesh used to construct the physics on |
| [in] | cycle | The simulation cycle (i.e. timestep iteration) to intialize the physics module to |
| [in] | time | The simulation time to initialize the physics module to |
| [in] | checkpoint_to_disk | A flag to save the transient states on disk instead of memory for the transient adjoint solves |
Definition at line 22 of file base_physics.cpp.
|
default |
Construct a new Base Physics object (copy constructor)
| other | The other base physics to copy from |
|
pure virtual |
Accessor for getting named finite element state adjoint solution from the physics modules.
| adjoint_name | The name of the Finite Element State adjoint solution to retrieve |
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
inlinevirtual |
Get a vector of the finite element state adjoint solution names.
Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 198 of file base_physics.hpp.
|
pure virtual |
Advance the state variables according to the chosen time integrator.
Advance the underlying ODE with the requested time integration scheme using the previously set timestep.
| dt | The increment of simulation time to advance the underlying physical system |
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
pure virtual |
Complete the setup and allocate the necessary data structures.
This finializes the underlying data structures in a solver and enables it to be run through a timestepping loop.
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanicsContact< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
inlinevirtual |
Compute the implicit sensitivity of the quantity of interest with respect to the initial condition fields.
reverseAdjointTimestep with an appropriate adjoint load must be called prior to this method as many times as the forward advance is called. Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 351 of file base_physics.hpp.
|
virtual |
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).
| parameter_index | the index of the parameter |
Reimplemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 105 of file base_physics.cpp.
|
virtual |
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).
Reimplemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanicsContact< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 111 of file base_physics.cpp.
|
virtual |
Get the current forward-solution cycle iteration number.
Definition at line 39 of file base_physics.cpp.
|
inlinevirtual |
Accessor for getting named finite element state dual (reaction) solution from the physics modules.
| dual_name | The name of the Finite Element State dual solution to retrieve |
Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 213 of file base_physics.hpp.
|
inlinevirtual |
Accessor for getting named finite element state dual adjoint (reaction adjoint load) from the physics modules.
| dual_name | The name of the Finite Element State dual (reaction adjoint load) solution to retrieve |
Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 227 of file base_physics.hpp.
|
inlinevirtual |
Get a vector of the finite element state dual (reaction) solution names.
Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 205 of file base_physics.hpp.
|
protected |
Accessor for getting all of the primal solutions from the physics modules at a given checkpointed cycle index.
| cycle | The cycle to retrieve state from |
Definition at line 374 of file base_physics.cpp.
|
virtual |
Get a timestep increment which has been previously checkpointed at the give cycle.
| cycle | The previous 'timestep' number where the timestep increment is requested |
Definition at line 393 of file base_physics.cpp.
|
protected |
Protected, non-virtual method to reset physics states to zero. This does not reset design parameters or shape.
| [in] | cycle | The simulation cycle (i.e. timestep iteration) to intialize the physics module to |
| [in] | time | The simulation time to initialize the physics module to |
Definition at line 57 of file base_physics.cpp.
|
virtual |
Initializes the Sidre structure for simulation summary data.
| [in] | datastore | Sidre DataStore where data are saved |
| [in] | t_final | Final time of the simulation |
| [in] | dt | The time step |
Definition at line 223 of file base_physics.cpp.
|
inline |
Check if the physics is setup as quasistatic.
Definition at line 129 of file base_physics.hpp.
|
inlinevirtual |
Accessor for getting a single named finite element dual solution from the physics modules at a given checkpointed cycle index.
| state_name | The name of the state to retrieve (e.g. "reaction") |
| cycle | The cycle to retrieve state from |
Definition at line 447 of file base_physics.hpp.
|
virtual |
Accessor for getting a single named finite element state primal solution from the physics modules at a given checkpointed cycle index.
| state_name | The name of the state to retrieve (e.g. "temperature", "displacement") |
| cycle | The cycle to retrieve state from |
Definition at line 340 of file base_physics.cpp.
| void serac::BasePhysics::loadCheckpointedStatesFromDisk | ( | int | cycle | ) |
load checkpointed states from disk into states array
| cycle | The cycle to retrieve state from |
Definition at line 365 of file base_physics.cpp.
|
virtual |
The maximum cycle (timestep iteration number) reached by the forward solver.
Definition at line 43 of file base_physics.cpp.
|
virtual |
Get the maximum time reached by the forward solver.
Definition at line 41 of file base_physics.cpp.
|
virtual |
Get the initial cycle (timestep iteration number) used by the forward solver.
Definition at line 47 of file base_physics.cpp.
|
virtual |
Get the initial time used by the forward solver.
Definition at line 45 of file base_physics.cpp.
|
virtual |
Output the current state of the PDE fields in Sidre format and optionally in Paraview format if paraview_output_dir is given.
| [in] | paraview_output_dir | Optional output directory for paraview visualization files |
Definition at line 187 of file base_physics.cpp.
|
inlinevirtual |
Accessor for getting named finite element state parameter fields from the physics modules.
| parameter_name | The name of the Finite Element State parameter to retrieve |
Definition at line 250 of file base_physics.hpp.
|
inlinevirtual |
Accessor for getting indexed finite element state parameter fields from the physics modules.
| parameter_index | The index of the Finite Element State parameter to retrieve |
Definition at line 272 of file base_physics.hpp.
|
inlinevirtual |
Get a vector of the finite element state parameter names.
Definition at line 287 of file base_physics.hpp.
|
pure virtual |
Base method to reset physics states to the initial time. This does not reset design parameters or shape.
| [in] | cycle | The simulation cycle (i.e. timestep iteration) to intialize the physics module to |
| [in] | time | The simulation time to initialize the physics module to |
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanicsContact< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
inlinevirtual |
Solve the adjoint reverse timestep problem.
Reimplemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 399 of file base_physics.hpp.
|
virtual |
Saves the summary data to the Sidre Datastore.
| [in] | datastore | Sidre DataStore where curves are saved |
| [in] | t | The current time of the simulation |
Definition at line 276 of file base_physics.cpp.
|
inlinevirtual |
Set the loads for the adjoint reverse timestep solve.
| string_to_dual | An unorder map from the state field name string to the finite element adjoint/dual load The adjoint load is d(qoi)/d(state) |
Reimplemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 372 of file base_physics.hpp.
|
inlinevirtual |
Set the dual loads (dirichlet values) for the adjoint reverse timestep solve This must be called after setAdjointLoad.
| string_to_bc | An unorder map from dual name string to finite element adjoint/state boundary condition The adjoint bc is d(qoi)/d(dual) |
Reimplemented in serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
Definition at line 386 of file base_physics.hpp.
|
virtual |
Deep copy a parameter field into the internally-owned parameter used for simulations.
| parameter_index | the index of the parameter |
| parameter_state | the values to use for the specified parameter |
The physics module constructs its own parameter FiniteElementState in the physics module constructor. This call sets the internally-owned parameter object by value (i.e. deep copies) from the given argument.
Definition at line 76 of file base_physics.cpp.
|
virtual |
Set the current shape displacement for the underlying mesh.
| shape_displacement | The shape displacement to copy for use in the physics module |
This updates the shape displacement field associated with the underlying mesh. Note that the input FiniteElementState is deep copied into the shape displacement object owned by the StateManager.
Definition at line 98 of file base_physics.cpp.
|
virtual |
Accessor for getting the shape displacement field from the physics modules.
Definition at line 96 of file base_physics.cpp.
|
protected |
Internally used accessor for getting the shape displacement sensitivity from the physics modules.
Definition at line 103 of file base_physics.cpp.
|
pure virtual |
Accessor for getting named finite element state primal solution from the physics modules.
| state_name | The name of the Finite Element State primal solution to retrieve |
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
pure virtual |
Get a vector of the finite element state primal solution names.
Implemented in serac::ThermomechanicsMonolithic< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, serac::Thermomechanics< order, dim, parameter_space >, serac::SolidMechanics< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >, and serac::HeatTransfer< order, dim, Parameters< parameter_space... >, std::integer_sequence< int, parameter_indices... > >.
|
virtual |
Get the current forward-solution time.
Definition at line 37 of file base_physics.cpp.
|
virtual |
Get a vector of the timestep sizes (i.e. \(\Delta t\)s) taken by the forward solver.
Definition at line 49 of file base_physics.cpp.
|
protected |
Update the paraview states, duals, parameters, and metadata (cycle, time) in preparation for output.
| paraview_output_dir | The directory to write the paraview output |
Definition at line 164 of file base_physics.cpp.
|
mutableprotected |
A container relating a checkpointed cycle and the associated finite element state fields.
Definition at line 621 of file base_physics.hpp.