29 namespace heat_transfer {
36 .relative_tol = 1.0e-6,
37 .absolute_tol = 1.0e-12,
38 .max_iterations = 200};
41 #ifdef MFEM_USE_STRUMPACK
51 .relative_tol = 1.0e-4,
52 .absolute_tol = 1.0e-8,
53 .max_iterations = 500,
81 template <
int order,
int dim,
typename parameters = Parameters<>,
82 typename parameter_indices = std::make_
integer_sequence<
int, parameters::n>>
86 template <
int order,
int dim,
typename... parameter_space,
int... parameter_indices>
91 static constexpr
int VALUE = 0, DERIVATIVE = 1;
92 static constexpr
int SHAPE = 0;
93 static constexpr
auto I = Identity<dim>();
98 static constexpr
auto NUM_STATE_VARS = 2;
119 std::shared_ptr<serac::Mesh> serac_mesh, std::vector<std::string> parameter_names = {},
int cycle = 0,
120 double time = 0.0,
bool checkpoint_to_disk =
false)
121 :
HeatTransfer(std::make_unique<EquationSolver>(nonlinear_opts, lin_opts, serac_mesh->getComm()),
122 timestepping_opts, physics_name, serac_mesh, parameter_names, cycle, time, checkpoint_to_disk)
143 const std::string& physics_name, std::shared_ptr<serac::Mesh> serac_mesh,
144 std::vector<std::string> parameter_names = {},
int cycle = 0,
double time = 0.0,
145 bool checkpoint_to_disk =
false)
146 :
BasePhysics(physics_name, serac_mesh, cycle, time, checkpoint_to_disk),
150 adjoint_temperature_(
151 StateManager::newState(H1<order>{}, detail::addPrefix(physics_name,
"adjoint_temperature"), mesh_->tag())),
152 implicit_sensitivity_temperature_start_of_step_(adjoint_temperature_.space(),
153 detail::addPrefix(physics_name,
"total_deriv_wrt_temperature")),
154 temperature_adjoint_load_(temperature_.space(), detail::addPrefix(physics_name,
"temperature_adjoint_load")),
155 temperature_rate_adjoint_load_(temperature_.space(),
156 detail::addPrefix(physics_name,
"temperature_rate_adjoint_load")),
157 residual_with_bcs_(temperature_.space().TrueVSize()),
158 nonlin_solver_(std::move(solver)),
159 ode_(temperature_.space().TrueVSize(),
160 {.time = time_, .u = u_, .dt = dt_, .du_dt = temperature_rate_, .previous_dt = previous_dt_},
161 *nonlin_solver_, bcs_)
164 mfemParMesh().Dimension() != dim,
165 axom::fmt::format(
"Compile time class dimension template parameter and runtime mesh dimension do not match"));
169 "EquationSolver argument is nullptr in HeatTransfer constructor. It is possible that it was previously moved.");
173 ode_.SetTimestepper(timestepping_opts.
timestepper);
175 is_quasistatic_ =
false;
177 is_quasistatic_ =
true;
180 states_.push_back(&temperature_);
181 if (!is_quasistatic_) {
182 states_.push_back(&temperature_rate_);
185 adjoints_.push_back(&adjoint_temperature_);
188 const mfem::ParFiniteElementSpace* test_space = &temperature_.space();
189 const mfem::ParFiniteElementSpace* shape_space = &mesh_->shapeDisplacementSpace();
191 std::array<
const mfem::ParFiniteElementSpace*,
sizeof...(parameter_space) + NUM_STATE_VARS> trial_spaces;
192 trial_spaces[0] = &temperature_.space();
193 trial_spaces[1] = &temperature_.space();
196 sizeof...(parameter_space) != parameter_names.size(),
197 axom::fmt::format(
"{} parameter spaces given in the template argument but {} parameter names were supplied.",
198 sizeof...(parameter_space), parameter_names.size()));
200 if constexpr (
sizeof...(parameter_space) > 0) {
201 tuple<parameter_space...> types{};
203 parameters_.emplace_back(mfemParMesh(), get<i>(types), detail::addPrefix(name_, parameter_names[i]));
205 trial_spaces[i + NUM_STATE_VARS] = &(parameters_[i].state->space());
210 std::make_unique<ShapeAwareFunctional<shape_trial, test(scalar_trial, scalar_trial, parameter_space...)>>(
211 shape_space, test_space, trial_spaces);
213 nonlin_solver_->setOperator(residual_with_bcs_);
215 int true_size = temperature_.space().TrueVSize();
216 u_.SetSize(true_size);
217 u_predicted_.SetSize(true_size);
219 initializeThermalStates();
232 std::shared_ptr<serac::Mesh> serac_mesh,
int cycle = 0,
double time = 0.0)
233 :
HeatTransfer(input_options.nonlin_solver_options, input_options.lin_solver_options,
234 input_options.timestepping_options, physics_name, serac_mesh, {}, cycle, time)
236 for (
const auto& mat : input_options.materials) {
237 if (std::holds_alternative<serac::heat_transfer::LinearIsotropicConductor>(mat)) {
238 setMaterial(std::get<serac::heat_transfer::LinearIsotropicConductor>(mat));
244 if (input_options.initial_temperature) {
245 auto temp = input_options.initial_temperature->constructScalar();
246 temperature_.project(*temp);
249 if (input_options.source_coef) {
252 SLIC_ERROR(
"'source' is not implemented yet in input files.");
256 std::map<std::string, input::BoundaryConditionInputOptions> sorted_bcs(input_options.boundary_conditions.begin(),
257 input_options.boundary_conditions.end());
258 for (
const auto& [bc_name, bc] : sorted_bcs) {
260 if (bc_name.find(
"temperature") != std::string::npos) {
261 std::shared_ptr<mfem::Coefficient> temp_coef(bc.coef_opts.constructScalar());
262 bcs_.addEssential(bc.attrs, temp_coef, temperature_.space(), *bc.coef_opts.component);
263 }
else if (bc_name.find(
"flux") != std::string::npos) {
266 SLIC_ERROR(
"'flux' is not implemented yet in input files.");
268 SLIC_WARNING_ROOT(
"Ignoring boundary condition with unknown name: " << physics_name);
280 time_end_step_ = 0.0;
284 temperature_rate_ = 0.0;
285 adjoint_temperature_ = 0.0;
286 implicit_sensitivity_temperature_start_of_step_ = 0.0;
287 temperature_adjoint_load_ = 0.0;
288 temperature_rate_adjoint_load_ = 0.0;
300 initializeThermalStates();
302 if (!checkpoint_to_disk_) {
303 checkpoint_states_.clear();
304 auto state_names = stateNames();
305 for (
const auto& state_name : state_names) {
306 checkpoint_states_[state_name].push_back(state(state_name));
319 void setTemperatureBCs(
const std::set<int>& temp_bdr, std::function<
double(
const mfem::Vector& x,
double t)> temp)
322 temp_bdr_coef_ = std::make_shared<mfem::FunctionCoefficient>(temp);
324 bcs_.addEssential(temp_bdr, temp_bdr_coef_, temperature_.space());
336 if (is_quasistatic_) {
340 for (
auto& bc : bcs_.essentials()) {
341 bc.setDofs(temperature_, time_);
343 nonlin_solver_->solve(temperature_);
354 double time_tmp = time_;
355 ode_.Step(temperature_, time_tmp, dt);
360 if (checkpoint_to_disk_) {
363 auto state_names = stateNames();
364 for (
const auto& state_name : state_names) {
365 checkpoint_states_[state_name].push_back(state(state_name));
369 if (cycle_ > max_cycle_) {
370 timesteps_.push_back(dt);
380 template <
typename MaterialType>
381 struct ThermalMaterialIntegrand {
392 template <
typename X,
typename T,
typename dT_dt,
typename... Params>
393 auto operator()(
double , X x, T temperature, dT_dt dtemp_dt, Params... params)
const
396 auto [u, du_dX] = temperature;
397 auto du_dt = get<VALUE>(dtemp_dt);
399 auto [heat_capacity, heat_flux] = material_(x, u, du_dX, params...);
401 return serac::tuple{heat_capacity * du_dt, -1.0 * heat_flux};
405 MaterialType material_;
432 template <
int... active_parameters,
typename MaterialType>
436 ThermalMaterialIntegrand<MaterialType>(material), domain);
440 template <
typename MaterialType>
456 mfem::FunctionCoefficient temp_coef(temp);
458 temp_coef.SetTime(time_);
459 temperature_.project(temp_coef);
487 template <
int... active_parameters,
typename SourceType>
490 residual_->AddDomainIntegral(
492 [source_function](
double t,
auto x,
auto temperature,
auto ,
auto... params) {
494 auto [u, du_dX] = temperature;
496 auto source = source_function(x, t, u, du_dX, params...);
505 template <
typename SourceType>
533 template <
int... active_parameters,
typename FluxType>
536 residual_->AddBoundaryIntegral(
538 [flux_function](
double t,
auto X,
auto u,
auto ,
auto... params) {
539 auto temp = get<VALUE>(u);
540 auto n =
cross(get<DERIVATIVE>(X));
542 return flux_function(X,
normalize(n), t, temp, params...);
548 template <
typename FluxType>
576 if (state_name ==
"temperature") {
578 }
else if (state_name ==
"temperature_rate") {
579 return temperature_rate_;
582 SLIC_ERROR_ROOT(axom::fmt::format(
"State '{}' requested from solid mechanics module '{}', but it doesn't exist",
598 if (state_name ==
"temperature") {
599 temperature_ = state;
600 if (!checkpoint_to_disk_) {
601 checkpoint_states_[
"temperature"][
static_cast<size_t>(cycle_)] = temperature_;
606 SLIC_ERROR_ROOT(axom::fmt::format(
607 "setState for state named '{}' requested from heat transfer module '{}', but it doesn't exist", state_name,
618 if (is_quasistatic_) {
619 return std::vector<std::string>{
"temperature"};
621 return std::vector<std::string>{
"temperature",
"temperature_rate"};
648 template <
int... active_parameters,
typename callable>
663 if (state_name ==
"temperature") {
664 return adjoint_temperature_;
667 SLIC_ERROR_ROOT(axom::fmt::format(
"Adjoint '{}' requested from solid mechanics module '{}', but it doesn't exist",
669 return adjoint_temperature_;
679 if (is_quasistatic_) {
681 temperature_.space().TrueVSize(),
683 [
this](
const mfem::Vector& u, mfem::Vector& r) {
684 const mfem::Vector res = (*residual_)(time_, shapeDisplacement(), u, temperature_rate_,
685 *parameters_[parameter_indices].state...);
691 r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
694 [
this](
const mfem::Vector& u) -> mfem::Operator& {
695 auto [r, drdu] = (*residual_)(time_, shapeDisplacement(), differentiate_wrt(u), temperature_rate_,
696 *parameters_[parameter_indices].state...);
698 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
703 temperature_.space().TrueVSize(),
705 [
this](
const mfem::Vector& du_dt, mfem::Vector& r) {
706 add(1.0, u_, dt_, du_dt, u_predicted_);
707 const mfem::Vector res =
708 (*residual_)(time_, shapeDisplacement(), u_predicted_, du_dt, *parameters_[parameter_indices].state...);
714 r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
717 [
this](
const mfem::Vector& du_dt) -> mfem::Operator& {
718 add(1.0, u_, dt_, du_dt, u_predicted_);
721 auto K = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), differentiate_wrt(u_predicted_),
722 du_dt, *parameters_[parameter_indices].state...));
723 std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
727 serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), u_predicted_, differentiate_wrt(du_dt),
728 *parameters_[parameter_indices].state...));
729 std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
732 J_.reset(mfem::Add(1.0, *m_mat, dt_, *k_mat));
733 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
739 if (checkpoint_to_disk_) {
742 checkpoint_states_.clear();
743 auto state_names = stateNames();
744 for (
const auto& state_name : state_names) {
745 checkpoint_states_[state_name].push_back(state(state_name));
764 virtual void setAdjointLoad(std::unordered_map<std::string, const serac::FiniteElementDual&> loads)
override
766 SLIC_ERROR_ROOT_IF(loads.size() == 0,
767 "Adjoint load container size must be greater than 0 in the heat transfer module.");
769 auto temp_adjoint_load = loads.find(
"temperature");
770 auto temp_rate_adjoint_load = loads.find(
"temperature_rate");
772 SLIC_ERROR_ROOT_IF(temp_adjoint_load == loads.end(),
"Adjoint load for \"temperature\" not found.");
774 temperature_adjoint_load_ = temp_adjoint_load->second;
776 temperature_adjoint_load_ *= -1.0;
778 if (temp_rate_adjoint_load != loads.end()) {
779 temperature_rate_adjoint_load_ = temp_rate_adjoint_load->second;
780 temperature_rate_adjoint_load_ *= -1.0;
816 auto& lin_solver = nonlin_solver_->linearSolver();
818 mfem::HypreParVector adjoint_essential(temperature_adjoint_load_);
819 adjoint_essential = 0.0;
821 SLIC_ERROR_ROOT_IF(cycle_ <= min_cycle_,
822 "Maximum number of adjoint timesteps exceeded! The number of adjoint timesteps must equal the "
823 "number of forward timesteps");
827 double dt = getCheckpointedTimestep(cycle_);
828 auto end_step_solution = getCheckpointedStates(cycle_ + 1);
830 temperature_ = end_step_solution.at(
"temperature");
832 if (is_quasistatic_) {
836 auto [_, drdu] = (*residual_)(time_, shapeDisplacement(),
differentiate_wrt(temperature_), temperature_rate_,
837 *parameters_[parameter_indices].state...);
838 auto jacobian = assemble(drdu);
839 auto J_T = std::unique_ptr<mfem::HypreParMatrix>(jacobian->Transpose());
841 for (
const auto& bc : bcs_.essentials()) {
842 bc.apply(*J_T, temperature_adjoint_load_, adjoint_essential);
845 lin_solver.SetOperator(*J_T);
846 lin_solver.Mult(temperature_adjoint_load_, adjoint_temperature_);
849 "Only backward Euler implemented for transient adjoint heat conduction.");
851 temperature_rate_ = end_step_solution.at(
"temperature_rate");
854 auto K = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(),
differentiate_wrt(temperature_),
855 temperature_rate_, *parameters_[parameter_indices].state...));
856 std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
859 auto M = serac::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), temperature_,
861 *parameters_[parameter_indices].state...));
862 std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
864 J_.reset(mfem::Add(1.0, *m_mat, dt, *k_mat));
865 auto J_T = std::unique_ptr<mfem::HypreParMatrix>(J_->Transpose());
869 mfem::HypreParVector modified_RHS(temperature_adjoint_load_);
871 modified_RHS.Add(1.0, temperature_rate_adjoint_load_);
872 modified_RHS.Add(-dt, implicit_sensitivity_temperature_start_of_step_);
874 for (
const auto& bc : bcs_.essentials()) {
875 bc.apply(*J_T, modified_RHS, adjoint_essential);
878 lin_solver.SetOperator(*J_T);
879 lin_solver.Mult(modified_RHS, adjoint_temperature_);
885 m_mat->Mult(adjoint_temperature_, implicit_sensitivity_temperature_start_of_step_);
886 implicit_sensitivity_temperature_start_of_step_ *= -1.0 / dt;
887 implicit_sensitivity_temperature_start_of_step_.Add(1.0 / dt,
888 temperature_rate_adjoint_load_);
891 time_end_step_ = time_;
908 auto drdparam = serac::get<DERIVATIVE>(d_residual_d_[parameter_field](time_end_step_));
909 auto drdparam_mat = assemble(drdparam);
911 drdparam_mat->MultTranspose(adjoint_temperature_, *parameters_[parameter_field].sensitivity);
913 return *parameters_[parameter_field].sensitivity;
927 serac::get<DERIVATIVE>((*residual_)(time_end_step_,
differentiate_wrt(shapeDisplacement()), temperature_,
928 temperature_rate_, *parameters_[parameter_indices].state...));
930 auto drdshape_mat = assemble(drdshape);
932 drdshape_mat->MultTranspose(adjoint_temperature_, shape_displacement_dual_);
934 return shapeDisplacementSensitivity();
947 return {{
"temperature", implicit_sensitivity_temperature_start_of_step_}};
985 std::unique_ptr<mfem::HypreParMatrix>
M_;
1007 std::unique_ptr<mfem::HypreParMatrix>
J_;
1011 std::unique_ptr<mfem::HypreParMatrix>
J_e_;
1032 std::array<std::function<decltype((*residual_)(
DifferentiateWRT<1>{}, 0.0, shape_displacement_, temperature_,
1033 temperature_rate_, *parameters_[parameter_indices].state...))(double)>,
1034 sizeof...(parameter_indices)>
1035 d_residual_d_ = {[&](
double TIME) {
1037 temperature_, temperature_rate_, *parameters_[parameter_indices].state...);
The base interface class for a generic PDE solver.
This is the abstract base class for a generic forward solver.
void initializeBasePhysicsStates(int cycle, double time)
Protected, non-virtual method to reset physics states to zero. This does not reset design parameters ...
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.
mfem_ext::StdFunctionOperator residual_with_bcs_
mfem::Operator that describes the weight residual and its gradient with respect to temperature
void setTemperature(const FiniteElementState temp)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void setSource(SourceType source_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
FiniteElementDual computeTimestepSensitivity(size_t parameter_field) override
Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoin...
const std::unordered_map< std::string, const serac::FiniteElementDual & > computeInitialConditionSensitivity() const override
Compute the implicit sensitivity of the quantity of interest with respect to the initial temperature.
void setFluxBCs(DependsOn< active_parameters... >, FluxType flux_function, Domain &domain)
Set the thermal flux boundary condition.
void addCustomBoundaryIntegral(DependsOn< active_parameters... >, callable qfunction, Domain &domain)
register a custom boundary integral calculation as part of the residual
const FiniteElementState & state(const std::string &state_name) const override
Accessor for getting named finite element state primal solution from the physics modules.
HeatTransfer(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< serac::Mesh > serac_mesh, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false)
Construct a new heat transfer object.
serac::FiniteElementDual temperature_rate_adjoint_load_
The downstream derivative of the qoi with repect to the primal temperature rate variable.
HeatTransfer(std::unique_ptr< serac::EquationSolver > solver, const serac::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< serac::Mesh > serac_mesh, std::vector< std::string > parameter_names={}, int cycle=0, double time=0.0, bool checkpoint_to_disk=false)
Construct a new heat transfer object.
double dt_
The current timestep.
void setMaterial(const MaterialType &material, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
virtual ~HeatTransfer()=default
Destroy the Thermal Solver object.
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.
serac::FiniteElementState adjoint_temperature_
The adjoint temperature finite element states, the multiplier on the residual for a given timestep.
virtual void setAdjointLoad(std::unordered_map< std::string, const serac::FiniteElementDual & > loads) override
Set the loads for the adjoint reverse timestep solve.
mfem_ext::FirstOrderODE ode_
the ordinary differential equation that describes how to solve for the time derivative of temperature...
void setTemperature(std::function< double(const mfem::Vector &x, double t)> temp)
Set the underlying finite element state to a prescribed temperature.
void setFluxBCs(FluxType flux_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
std::unique_ptr< EquationSolver > nonlin_solver_
the specific methods and tolerances specified to solve the nonlinear residual equations
mfem::Vector u_
Predicted temperature true dofs.
const serac::FiniteElementState & temperature() const
Get the temperature state.
void setSource(DependsOn< active_parameters... >, SourceType source_function, Domain &domain)
Set the thermal source function.
mfem::Vector u_predicted_
Predicted temperature true dofs.
const FiniteElementDual & computeTimestepShapeSensitivity() override
Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoin...
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)
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 > J_e_
void initializeThermalStates()
Non virtual method to reset thermal states to zero. This does not reset design parameters or shape.
virtual std::vector< std::string > stateNames() const override
Get a vector of the finite element state primal solution names.
const serac::FiniteElementState & temperatureRate() const
Get the temperature rate of change state.
void reverseAdjointTimestep() override
Solve the adjoint problem.
serac::FiniteElementDual temperature_adjoint_load_
The downstream derivative of the qoi with repect to the primal temperature variable.
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 temper...
HeatTransfer(const HeatTransferInputOptions &input_options, const std::string &physics_name, std::shared_ptr< serac::Mesh > serac_mesh, int cycle=0, double time=0.0)
Construct a new Nonlinear HeatTransfer Solver object.
void setMaterial(DependsOn< active_parameters... >, const MaterialType &material, Domain &domain)
Set the thermal material model for the physics solver.
serac::FiniteElementState temperature_
The temperature finite element state.
void completeSetup() override
Complete the initialization and allocation of the data structures.
std::unique_ptr< mfem::HypreParMatrix > J_
Assembled sparse matrix for the Jacobian.
void advanceTimestep(double dt) override
Advance the heat conduction physics module in time.
const FiniteElementState & adjoint(const std::string &state_name) const override
Accessor for getting named finite element state adjoint solution from the physics modules.
FiniteElementState temperature_rate_
Rate of change in temperature at the current adjoint timestep.
double time_end_step_
End of step time used in reverse mode so that the time can be decremented on reverse steps.
double previous_dt_
The previous timestep.
void setState(const std::string &state_name, const FiniteElementState &state) override
Set the primal solution field (temperature) for the underlying heat transfer solver.
std::shared_ptr< mfem::Coefficient > temp_bdr_coef_
Coefficient containing the essential boundary values.
std::unique_ptr< mfem::HypreParMatrix > M_
Assembled mass matrix.
An object containing the solver for a heat transfer PDE.
static FiniteElementState newState(FunctionSpace space, const std::string &state_name, const std::string &mesh_tag)
Factory method for creating a new FEState object.
FirstOrderODE is a class wrapping mfem::TimeDependentOperator so that the user can use std::function ...
StdFunctionOperator is a class wrapping mfem::Operator so that the user can use std::function to defi...
A file defining some enums and structs that are used by the different physics modules.
Serac mesh class which assists in constructing the appropriate parallel mfem meshes and registering a...
const TimesteppingOptions default_timestepping_options
Reasonable defaults for dynamic heat transfer simulations.
const TimesteppingOptions default_static_options
Reasonable defaults for static heat transfer simulations.
const LinearSolverOptions default_linear_options
Reasonable defaults for most thermal linear solver options.
const NonlinearSolverOptions default_nonlinear_options
Reasonable defaults for most thermal nonlinear solver options.
const LinearSolverOptions direct_linear_options
the default direct solver option for solving the linear stiffness equations
Accelerator functionality.
auto differentiate_wrt(const mfem::Vector &v)
this function is intended to only be used in combination with serac::Functional::operator(),...
tuple(T...) -> tuple< T... >
Class template argument deduction rule for tuples.
auto cross(const tensor< T, 3, 2 > &A)
compute the cross product of the columns of A: A(:,1) x A(:,2)
SERAC_HOST_DEVICE auto normalize(const tensor< T, n... > &A)
Normalizes the tensor Each element is divided by the Frobenius norm of the tensor,...
bool holds_alternative(const variant< T0, T1 > &v)
Checks whether a variant's active member is of a certain type.
Wrapper classes for using MFEM's ODE solver objects.
Wrapper of serac::Functional for evaluating integrals and derivatives of quantities with shape displa...
This file contains the declaration of the StateManager class.
Classes for defining an mfem::Operator with std::functions.
Compile-time alias for a dimension.
a class for representing a geometric region that can be used for integration
auto operator()(double, X x, T temperature, dT_dt dtemp_dt, Params... params) const
Evaluate integrand.
ThermalMaterialIntegrand(MaterialType material)
Construct a ThermalMaterialIntegrand functor with material model of type MaterialType.
Parameters for an iterative linear solution scheme.
LinearSolver linear_solver
Linear solver selection.
Nonlinear solution scheme parameters.
NonlinearSolver nonlin_solver
Nonlinear solver selection.
a struct that is used in the physics modules to clarify which template arguments are user-controlled ...
A timestep and boundary condition enforcement method for a dynamic solver.
TimestepMethod timestepper
The timestepping method to be applied.
DirichletEnforcementMethod enforcement_method
The essential boundary enforcement method to use.
Linear anisotropic thermal material model.
This is a class that mimics most of std::tuple's interface, except that it is usable in CUDA kernels ...
A sentinel struct for eliding no-op tensor operations.
The material and load types for the thermal functional physics module.