20 #include <unordered_map>
25 #include "axom/sidre.hpp"
26 #include "axom/fmt.hpp"
27 #include "axom/slic.hpp"
47 std::string addPrefix(
const std::string& prefix,
const std::string& target);
54 std::string removePrefix(
const std::string& prefix,
const std::string& target);
73 bool checkpoint_to_disk =
false);
87 virtual double time()
const;
94 virtual int cycle()
const;
101 virtual double maxTime()
const;
108 virtual double minTime()
const;
136 virtual const std::vector<double>&
timesteps()
const;
205 virtual std::vector<std::string>
dualNames()
const {
return {}; }
215 SLIC_ERROR_ROOT(axom::fmt::format(
"dual '{}' requested from physics module '{}' which does not support duals",
229 SLIC_ERROR_ROOT(axom::fmt::format(
230 "dualAdjoint '{}' requested from physics module '{}' which does not support duals", dual_name,
name_));
252 std::string appended_name = detail::addPrefix(
name_, parameter_name);
260 SLIC_ERROR_ROOT(axom::fmt::format(
"Parameter {} requested from physics module {}, but it doesn't exist.",
261 parameter_name,
name_));
276 axom::fmt::format(
"Parameter index {} requested, but only {} parameters exist in physics module {}.",
289 std::vector<std::string> parameter_names;
295 return parameter_names;
354 SLIC_WARNING_ROOT(axom::fmt::format(
"Initial condition sensitivities not enabled in physics module {}",
name_));
372 virtual void setAdjointLoad(std::unordered_map<std::string, const serac::FiniteElementDual&> string_to_dual)
374 if (!string_to_dual.empty()) {
376 axom::fmt::format(
"Failed to setAdjointLoad. Adjoint analysis not defined for physics module {}",
name_));
386 virtual void setDualAdjointBcs(std::unordered_map<std::string, const serac::FiniteElementState&> string_to_bc)
388 if (!string_to_bc.empty()) {
390 axom::fmt::format(
"Failed to setDualAdjointBCs. Adjoint analysis not defined for physics module {}",
name_));
401 SLIC_ERROR_ROOT(axom::fmt::format(
"Adjoint analysis not defined for physics module {}",
name_));
420 virtual void outputStateToDisk(std::optional<std::string> paraview_output_dir = {})
const;
448 [[maybe_unused]]
int cycle)
450 SLIC_ERROR_ROOT(axom::fmt::format(
"loadCheckpointedDual not enabled in physics module {}",
name_));
468 virtual void initializeSummary(axom::sidre::DataStore& datastore,
const double t_final,
const double dt)
const;
476 virtual void saveSummary(axom::sidre::DataStore& datastore,
const double t)
const;
557 std::vector<const serac::FiniteElementState*>
states_;
562 std::vector<const serac::FiniteElementState*>
adjoints_;
567 std::vector<const serac::FiniteElementDual*>
duals_;
584 template <
typename FunctionSpace>
587 state = std::make_unique<FiniteElementState>(
mesh, space,
name);
595 std::unique_ptr<serac::FiniteElementState>
state;
This file contains the declaration of the boundary condition manager class.
This is the abstract base class for a generic forward solver.
std::optional< int > cached_checkpoint_cycle_
An optional int for disk-based checkpointing containing the cycle number of the last retrieved checkp...
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 ...
virtual void completeSetup()=0
Complete the setup and allocate the necessary data structures.
std::unordered_map< std::string, serac::FiniteElementState > cached_checkpoint_states_
A container relating a checkpointed cycle and the associated finite element state fields.
virtual void setState(const std::string &, const FiniteElementState &)=0
Set the primal solution field values of the underlying physics solver.
std::vector< const serac::FiniteElementDual * > duals_
List of finite element duals associated with this physics module.
virtual void initializeSummary(axom::sidre::DataStore &datastore, const double t_final, const double dt) const
Initializes the Sidre structure for simulation summary data.
const mfem::ParMesh & mfemParMesh() const
Returns a reference to the mfem ParMesh object.
virtual void saveSummary(axom::sidre::DataStore &datastore, const double t) const
Saves the summary data to the Sidre Datastore.
double max_time_
The maximum time reached for the forward solver.
virtual std::vector< std::string > stateNames() const =0
Get a vector of the finite element state primal solution names.
int cycle_
Current cycle (forward pass time iteration count)
std::vector< const serac::FiniteElementState * > states_
List of finite element primal states associated with this physics module.
virtual void initializationStep()
Initialize any fields nessary for before the first step of the time integration.
serac::FiniteElementState shape_displacement_
The shape displacement field.
virtual double time() const
Get the current forward-solution time.
std::string name() const
Return the name of the physics.
virtual ~BasePhysics()=default
Destroy the Base Solver object.
std::unordered_map< std::string, std::vector< serac::FiniteElementState > > checkpoint_states_
A map containing optionally in-memory checkpointed primal states for transient adjoint solvers.
virtual FiniteElementDual computeTimestepSensitivity(size_t parameter_index)
Compute the implicit sensitivity of the quantity of interest used in defining the adjoint load with r...
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 afte...
double min_time_
The time the forward solver was initialized to.
virtual const FiniteElementState & parameter(const std::string ¶meter_name) const
Accessor for getting named finite element state parameter fields from the physics modules.
virtual std::vector< std::string > adjointNames() const
Get a vector of the finite element state adjoint solution names.
virtual std::vector< std::string > dualNames() const
Get a vector of the finite element state dual (reaction) solution names.
bool isQuasistatic() const
Check if the physics is setup as quasistatic.
virtual void reverseAdjointInitializationStep()
Compute adjoint sensitivities back through initializationStep.
serac::FiniteElementDual shape_displacement_dual_
The shape displacement field sensitivity.
double time_
Current time for the forward pass.
virtual void reverseAdjointTimestep()
Solve the adjoint reverse timestep problem.
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 sha...
virtual void resetAdjointStates()
Base method to reset physics states back to the end of time to start adjoint calculations again....
virtual const FiniteElementState & adjoint(const std::string &adjoint_name) const =0
Accessor for getting named finite element state adjoint solution from the physics modules.
int min_cycle_
The cycle the forward solver was initialized to.
bool checkpoint_to_disk_
A flag denoting whether to save the state to disk or memory as needed for dynamic adjoint solves.
double dt_
Current time step.
double ode_time_point_
The value of time at which the ODE solver wants to evaluate the residual.
int max_cycle_
The maximum cycle (forward pass iteration count) reached by the forward solver.
virtual int cycle() const
Get the current forward-solution cycle iteration number.
virtual void setShapeDisplacement(const FiniteElementState &shape_displacement)
Set the current shape displacement for the underlying mesh.
std::vector< const serac::FiniteElementState * > dual_adjoints_
List of adjoint finite element duals associated with this physics module.
void loadCheckpointedStatesFromDisk(int cycle)
load checkpointed states from disk into states array
virtual const FiniteElementState & shapeDisplacement() const
Accessor for getting the shape displacement field from the physics modules.
std::shared_ptr< serac::Mesh > mesh_
The primary mesh.
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 cyc...
BoundaryConditionManager bcs_
Boundary condition manager instance.
bool is_quasistatic_
Whether the simulation is time-independent.
BasePhysics(BasePhysics &&other)=default
Construct a new Base Physics object (copy constructor)
void UpdateParaviewDataCollection(const std::string ¶view_output_dir) const
Update the paraview states, duals, parameters, and metadata (cycle, time) in preparation for output.
virtual const std::vector< double > & timesteps() const
Get a vector of the timestep sizes (i.e. s) taken by the forward solver.
virtual const FiniteElementState & state(const std::string &state_name) const =0
Accessor for getting named finite element state primal solution from the physics modules.
virtual double minTime() const
Get the initial time used by the forward solver.
virtual void advanceTimestep(double dt)=0
Advance the state variables according to the chosen time integrator.
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 fi...
std::unique_ptr< mfem::ParaViewDataCollection > paraview_dc_
DataCollection pointer for optional 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.
virtual double getCheckpointedTimestep(int cycle) const
Get a timestep increment which has been previously checkpointed at the give cycle.
virtual const FiniteElementDual & dual(const std::string &dual_name) const
Accessor for getting named finite element state dual (reaction) solution from the physics modules.
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 parav...
virtual std::vector< std::string > parameterNames() const
Get a vector of the finite element state parameter names.
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.
std::vector< ParameterInfo > parameters_
A vector of the parameters associated with this physics module.
const serac::Mesh & mesh() const
Returns a reference to the mesh object.
virtual void setAdjointLoad(std::unordered_map< std::string, const serac::FiniteElementDual & > string_to_dual)
Set the loads for the adjoint reverse timestep solve.
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.
std::string name_
Name of the physics module.
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...
virtual const FiniteElementDual & computeTimestepShapeSensitivity()
Compute the implicit sensitivity of the quantity of interest used in defining the adjoint load with r...
virtual int minCycle() const
Get the initial cycle (timestep iteration number) used by the forward solver.
void CreateParaviewDataCollection() const
Create a paraview data collection for the physics package if requested.
virtual const FiniteElementState & parameter(std::size_t parameter_index) const
Accessor for getting indexed finite element state parameter fields from the physics modules.
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 ...
std::vector< const serac::FiniteElementState * > adjoints_
List of finite element adjoint states associated with this physics module.
virtual double maxTime() const
Get the maximum time reached by the forward solver.
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.
void initializeBasePhysicsStates(int cycle, double time)
Protected, non-virtual method to reset physics states to zero. This does not reset design parameters ...
MPI_Comm comm_
The MPI communicator.
std::vector< double > timesteps_
A vector of the timestep sizes (i.e. ) taken by the forward solver.
static constexpr int FLOAT_PRECISION_
Number of significant figures to output for floating-point.
const FiniteElementDual & shapeDisplacementSensitivity() const
Internally used accessor for getting the shape displacement sensitivity from the physics modules.
virtual int maxCycle() const
The maximum cycle (timestep iteration number) reached by the forward solver.
A container for the boundary condition information relating to a specific physics module.
Class for encapsulating the dual vector space of a finite element space (i.e. the space of linear for...
Class for encapsulating the critical MFEM components of a primal finite element field.
std::string name() const
Returns the name of the FEState (field)
Helper class for constructing a mesh consistent with serac.
static void storeDual(FiniteElementDual &dual)
Store a pre-constructed finite element dual in the state manager.
static void storeState(FiniteElementState &state)
Store a pre-constructed finite element state in the state manager.
A file defining some enums and structs that are used by the different physics modules.
This file contains the declaration of an equation solver wrapper.
This contains a class that represents the dual of a finite element vector space, i....
This file contains the declaration of structure that manages the MFEM objects that make up the state ...
Accelerator functionality.
This file contains the declaration of the StateManager class.
The information needed for the physics parameters stored as Finite Element State fields.
std::unique_ptr< serac::FiniteElementState > state
The finite element states representing user-defined and owned parameter fields.
std::unique_ptr< serac::FiniteElementState > previous_state
The finite element state representing the parameter at the previous evaluation.
std::unique_ptr< serac::FiniteElementDual > sensitivity
The sensitivities (dual vectors) of the QOI encoded in the adjoint load with respect to each of the i...
ParameterInfo(mfem::ParMesh &mesh, FunctionSpace space, const std::string &name="")
Construct a new Parameter Info object.
a small POD class for tracking function space metadata