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<smith::Mesh> smith_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, smith_mesh->getComm()),
122 timestepping_opts, physics_name, smith_mesh, parameter_names, cycle, time, checkpoint_to_disk)
143 const std::string& physics_name, std::shared_ptr<smith::Mesh> smith_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, smith_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 std::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 std::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<smith::Mesh> smith_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, smith_mesh, {}, cycle, time)
236 for (
const auto& mat : input_options.materials) {
237 if (std::holds_alternative<smith::heat_transfer::LinearIsotropicConductor>(mat)) {
238 setMaterial(std::get<smith::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 template <
typename AppliedTemperatureFunction>
322 auto mfem_coefficient_function = ([applied_temperature](
const mfem::Vector& X_mfem,
double t) {
323 auto X = make_tensor<dim>([&X_mfem](
int k) {
return X_mfem[k]; });
324 return applied_temperature(X, t);
326 temp_bdr_coef_ = std::make_shared<mfem::FunctionCoefficient>(mfem_coefficient_function);
327 auto dof_list = domain.
dof_list(&temperature_.space());
328 bcs_.addEssential(dof_list, temp_bdr_coef_, temperature_.space(), 0);
339 void setTemperatureBCs(
const std::set<int>& temp_bdr, std::function<
double(
const mfem::Vector& x,
double t)> temp)
342 temp_bdr_coef_ = std::make_shared<mfem::FunctionCoefficient>(temp);
344 bcs_.addEssential(temp_bdr, temp_bdr_coef_, temperature_.space());
356 if (is_quasistatic_) {
360 for (
auto& bc : bcs_.essentials()) {
361 bc.setDofs(temperature_, time_);
363 nonlin_solver_->solve(temperature_);
374 double time_tmp = time_;
375 ode_.Step(temperature_, time_tmp, dt);
380 if (checkpoint_to_disk_) {
383 auto state_names = stateNames();
384 for (
const auto& state_name : state_names) {
385 checkpoint_states_[state_name].push_back(state(state_name));
389 if (cycle_ > max_cycle_) {
390 timesteps_.push_back(dt);
400 template <
typename MaterialType>
401 struct ThermalMaterialIntegrand {
412 template <
typename X,
typename T,
typename dT_dt,
typename... Params>
413 auto operator()(
double , X x, T temperature, dT_dt dtemp_dt, Params... params)
const
416 auto [u, du_dX] = temperature;
417 auto du_dt = get<VALUE>(dtemp_dt);
419 auto [heat_capacity, heat_flux] = material_(x, u, du_dX, params...);
421 return smith::tuple{heat_capacity * du_dt, -1.0 * heat_flux};
425 MaterialType material_;
452 template <
int... active_parameters,
typename MaterialType>
456 ThermalMaterialIntegrand<MaterialType>(material), domain);
460 template <
typename MaterialType>
476 mfem::FunctionCoefficient temp_coef(temp);
478 temp_coef.SetTime(time_);
479 temperature_.project(temp_coef);
507 template <
int... active_parameters,
typename SourceType>
510 residual_->AddDomainIntegral(
512 [source_function](
double t,
auto x,
auto temperature,
auto ,
auto... params) {
514 auto [u, du_dX] = temperature;
516 auto source = source_function(x, t, u, du_dX, params...);
525 template <
typename SourceType>
553 template <
int... active_parameters,
typename FluxType>
556 residual_->AddBoundaryIntegral(
558 [flux_function](
double t,
auto X,
auto u,
auto ,
auto... params) {
559 auto temp = get<VALUE>(u);
560 auto n =
cross(get<DERIVATIVE>(X));
562 return flux_function(X,
normalize(n), t, temp, params...);
568 template <
typename FluxType>
596 if (state_name ==
"temperature") {
598 }
else if (state_name ==
"temperature_rate") {
599 return temperature_rate_;
603 std::format(
"State '{}' requested from solid mechanics module '{}', but it doesn't exist", state_name, name_));
618 if (state_name ==
"temperature") {
619 temperature_ = state;
620 if (!checkpoint_to_disk_) {
621 checkpoint_states_[
"temperature"][
static_cast<size_t>(cycle_)] = temperature_;
627 std::format(
"setState for state named '{}' requested from heat transfer module '{}', but it doesn't exist",
638 if (is_quasistatic_) {
639 return std::vector<std::string>{
"temperature"};
641 return std::vector<std::string>{
"temperature",
"temperature_rate"};
668 template <
int... active_parameters,
typename callable>
683 if (state_name ==
"temperature") {
684 return adjoint_temperature_;
687 SLIC_ERROR_ROOT(std::format(
"Adjoint '{}' requested from solid mechanics module '{}', but it doesn't exist",
689 return adjoint_temperature_;
699 if (is_quasistatic_) {
701 temperature_.space().TrueVSize(),
703 [
this](
const mfem::Vector& u, mfem::Vector& r) {
704 const mfem::Vector res = (*residual_)(time_, shapeDisplacement(), u, temperature_rate_,
705 *parameters_[parameter_indices].state...);
711 r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
714 [
this](
const mfem::Vector& u) -> mfem::Operator& {
715 auto [r, drdu] = (*residual_)(time_, shapeDisplacement(), differentiate_wrt(u), temperature_rate_,
716 *parameters_[parameter_indices].state...);
718 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
723 temperature_.space().TrueVSize(),
725 [
this](
const mfem::Vector& du_dt, mfem::Vector& r) {
726 add(1.0, u_, dt_, du_dt, u_predicted_);
727 const mfem::Vector res =
728 (*residual_)(time_, shapeDisplacement(), u_predicted_, du_dt, *parameters_[parameter_indices].state...);
734 r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
737 [
this](
const mfem::Vector& du_dt) -> mfem::Operator& {
738 add(1.0, u_, dt_, du_dt, u_predicted_);
741 auto K = smith::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), differentiate_wrt(u_predicted_),
742 du_dt, *parameters_[parameter_indices].state...));
743 std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
747 smith::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), u_predicted_, differentiate_wrt(du_dt),
748 *parameters_[parameter_indices].state...));
749 std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
752 J_.reset(mfem::Add(1.0, *m_mat, dt_, *k_mat));
753 J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
759 if (checkpoint_to_disk_) {
762 checkpoint_states_.clear();
763 auto state_names = stateNames();
764 for (
const auto& state_name : state_names) {
765 checkpoint_states_[state_name].push_back(state(state_name));
784 virtual void setAdjointLoad(std::unordered_map<std::string, const smith::FiniteElementDual&> loads)
override
786 SLIC_ERROR_ROOT_IF(loads.size() == 0,
787 "Adjoint load container size must be greater than 0 in the heat transfer module.");
789 auto temp_adjoint_load = loads.find(
"temperature");
790 auto temp_rate_adjoint_load = loads.find(
"temperature_rate");
792 SLIC_ERROR_ROOT_IF(temp_adjoint_load == loads.end(),
"Adjoint load for \"temperature\" not found.");
794 temperature_adjoint_load_ = temp_adjoint_load->second;
796 temperature_adjoint_load_ *= -1.0;
798 if (temp_rate_adjoint_load != loads.end()) {
799 temperature_rate_adjoint_load_ = temp_rate_adjoint_load->second;
800 temperature_rate_adjoint_load_ *= -1.0;
836 auto& lin_solver = nonlin_solver_->linearSolver();
838 mfem::HypreParVector adjoint_essential(temperature_adjoint_load_);
839 adjoint_essential = 0.0;
841 SLIC_ERROR_ROOT_IF(cycle_ <= min_cycle_,
842 "Maximum number of adjoint timesteps exceeded! The number of adjoint timesteps must equal the "
843 "number of forward timesteps");
847 double dt = getCheckpointedTimestep(cycle_);
848 auto end_step_solution = getCheckpointedStates(cycle_ + 1);
850 temperature_ = end_step_solution.at(
"temperature");
852 if (is_quasistatic_) {
856 auto [_, drdu] = (*residual_)(time_, shapeDisplacement(),
differentiate_wrt(temperature_), temperature_rate_,
857 *parameters_[parameter_indices].state...);
858 auto jacobian = assemble(drdu);
859 auto J_T = std::unique_ptr<mfem::HypreParMatrix>(jacobian->Transpose());
861 for (
const auto& bc : bcs_.essentials()) {
862 bc.apply(*J_T, temperature_adjoint_load_, adjoint_essential);
865 lin_solver.SetOperator(*J_T);
866 lin_solver.Mult(temperature_adjoint_load_, adjoint_temperature_);
869 "Only backward Euler implemented for transient adjoint heat conduction.");
871 temperature_rate_ = end_step_solution.at(
"temperature_rate");
874 auto K = smith::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(),
differentiate_wrt(temperature_),
875 temperature_rate_, *parameters_[parameter_indices].state...));
876 std::unique_ptr<mfem::HypreParMatrix> k_mat(assemble(K));
879 auto M = smith::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), temperature_,
881 *parameters_[parameter_indices].state...));
882 std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
884 J_.reset(mfem::Add(1.0, *m_mat, dt, *k_mat));
885 auto J_T = std::unique_ptr<mfem::HypreParMatrix>(J_->Transpose());
889 mfem::HypreParVector modified_RHS(temperature_adjoint_load_);
891 modified_RHS.Add(1.0, temperature_rate_adjoint_load_);
892 modified_RHS.Add(-dt, implicit_sensitivity_temperature_start_of_step_);
894 for (
const auto& bc : bcs_.essentials()) {
895 bc.apply(*J_T, modified_RHS, adjoint_essential);
898 lin_solver.SetOperator(*J_T);
899 lin_solver.Mult(modified_RHS, adjoint_temperature_);
905 m_mat->Mult(adjoint_temperature_, implicit_sensitivity_temperature_start_of_step_);
906 implicit_sensitivity_temperature_start_of_step_ *= -1.0 / dt;
907 implicit_sensitivity_temperature_start_of_step_.Add(1.0 / dt,
908 temperature_rate_adjoint_load_);
911 time_end_step_ = time_;
928 auto drdparam = smith::get<DERIVATIVE>(d_residual_d_[parameter_field](time_end_step_));
929 auto drdparam_mat = assemble(drdparam);
931 drdparam_mat->MultTranspose(adjoint_temperature_, *parameters_[parameter_field].sensitivity);
933 return *parameters_[parameter_field].sensitivity;
947 smith::get<DERIVATIVE>((*residual_)(time_end_step_,
differentiate_wrt(shapeDisplacement()), temperature_,
948 temperature_rate_, *parameters_[parameter_indices].state...));
950 auto drdshape_mat = assemble(drdshape);
952 drdshape_mat->MultTranspose(adjoint_temperature_, shape_displacement_dual_);
954 return shapeDisplacementSensitivity();
967 return {{
"temperature", implicit_sensitivity_temperature_start_of_step_}};
1005 std::unique_ptr<mfem::HypreParMatrix>
M_;
1027 std::unique_ptr<mfem::HypreParMatrix>
J_;
1031 std::unique_ptr<mfem::HypreParMatrix>
J_e_;
1052 std::array<std::function<decltype((*residual_)(
DifferentiateWRT<1>{}, 0.0, shape_displacement_, temperature_,
1053 temperature_rate_, *parameters_[parameter_indices].state...))(double)>,
1054 sizeof...(parameter_indices)>
1055 d_residual_d_ = {[&](
double TIME) {
1057 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.
void setMaterial(DependsOn< active_parameters... >, const MaterialType &material, Domain &domain)
Set the thermal material model for the physics solver.
void setTemperatureBCs(AppliedTemperatureFunction applied_temperature, Domain &domain)
Set essential temperature boundary conditions (strongly enforced)
void setFluxBCs(FluxType flux_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
smith::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...
const FiniteElementState & adjoint(const std::string &state_name) const override
Accessor for getting named finite element state adjoint solution from the physics modules.
void setSource(DependsOn< active_parameters... >, SourceType source_function, Domain &domain)
Set the thermal source function.
const FiniteElementDual & computeTimestepShapeSensitivity() override
Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoin...
double time_end_step_
End of step time used in reverse mode so that the time can be decremented on reverse steps.
const std::unordered_map< std::string, const smith::FiniteElementDual & > computeInitialConditionSensitivity() const override
Compute the implicit sensitivity of the quantity of interest with respect to the initial temperature.
smith::FiniteElementDual temperature_rate_adjoint_load_
The downstream derivative of the qoi with repect to the primal temperature rate variable.
const smith::FiniteElementState & temperatureRate() const
Get the temperature rate of change state.
FiniteElementDual computeTimestepSensitivity(size_t parameter_field) override
Compute the implicit sensitivity of the quantity of interest used in defining the load for the adjoin...
double previous_dt_
The previous timestep.
std::unique_ptr< ShapeAwareFunctional< shape_trial, test(scalar_trial, scalar_trial, parameter_space...)> > residual_
smith::Functional that is used to calculate the residual and its derivatives
void addCustomBoundaryIntegral(DependsOn< active_parameters... >, callable qfunction, Domain &domain)
register a custom boundary integral calculation as part of the residual
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)
void setState(const std::string &state_name, const FiniteElementState &state) override
Set the primal solution field (temperature) for the underlying heat transfer solver.
HeatTransfer(std::unique_ptr< smith::EquationSolver > solver, const smith::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_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.
mfem::Vector u_
Predicted temperature true dofs.
void advanceTimestep(double dt) override
Advance the heat conduction physics module in time.
mfem::Vector u_predicted_
Predicted temperature true dofs.
mfem_ext::FirstOrderODE ode_
the ordinary differential equation that describes how to solve for the time derivative of temperature...
smith::FiniteElementState temperature_
The temperature finite element state.
void initializeThermalStates()
Non virtual method to reset thermal states to zero. This does not reset design parameters or shape.
HeatTransfer(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts, const smith::TimesteppingOptions timestepping_opts, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_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.
std::unique_ptr< mfem::HypreParMatrix > J_
Assembled sparse matrix for the Jacobian.
std::shared_ptr< mfem::Coefficient > temp_bdr_coef_
Coefficient containing the essential boundary values.
void setMaterial(const MaterialType &material, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
virtual void setAdjointLoad(std::unordered_map< std::string, const smith::FiniteElementDual & > loads) override
Set the loads for the adjoint reverse timestep solve.
void setFluxBCs(DependsOn< active_parameters... >, FluxType flux_function, Domain &domain)
Set the thermal flux boundary condition.
mfem_ext::StdFunctionOperator residual_with_bcs_
mfem::Operator that describes the weight residual and its gradient with respect to temperature
const FiniteElementState & state(const std::string &state_name) const override
Accessor for getting named finite element state primal solution from the physics modules.
void completeSetup() override
Complete the initialization and allocation of the data structures.
smith::FiniteElementState adjoint_temperature_
The adjoint temperature finite element states, the multiplier on the residual for a given timestep.
double dt_
The current timestep.
void setSource(SourceType source_function, Domain &domain)
This is an overloaded member function, provided for convenience. It differs from the above function o...
HeatTransfer(const HeatTransferInputOptions &input_options, const std::string &physics_name, std::shared_ptr< smith::Mesh > smith_mesh, int cycle=0, double time=0.0)
Construct a new Nonlinear HeatTransfer Solver object.
void reverseAdjointTimestep() override
Solve the adjoint problem.
const smith::FiniteElementState & temperature() const
Get the temperature state.
std::unique_ptr< EquationSolver > nonlin_solver_
the specific methods and tolerances specified to solve the nonlinear residual equations
FiniteElementState temperature_rate_
Rate of change in temperature at the current adjoint timestep.
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.
void setTemperature(std::function< double(const mfem::Vector &x, double t)> temp)
Set the underlying finite element state to a prescribed temperature.
std::unique_ptr< mfem::HypreParMatrix > M_
Assembled mass matrix.
void setTemperature(const FiniteElementState temp)
This is an overloaded member function, provided for convenience. It differs from the above function o...
smith::FiniteElementDual temperature_adjoint_load_
The downstream derivative of the qoi with repect to the primal temperature variable.
virtual ~HeatTransfer()=default
Destroy the Thermal Solver object.
virtual std::vector< std::string > stateNames() const override
Get a vector of the finite element state primal solution names.
std::unique_ptr< mfem::HypreParMatrix > J_e_
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.
Smith 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 direct_linear_options
the default direct solver option for solving the linear stiffness equations
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.
Accelerator functionality.
tuple(T...) -> tuple< T... >
Class template argument deduction rule for tuples.
SMITH_HOST_DEVICE auto cross(const tensor< T, 3, 2 > &A)
compute the cross product of the columns of A: A(:,1) x A(:,2)
auto differentiate_wrt(const mfem::Vector &v)
this function is intended to only be used in combination with smith::Functional::operator(),...
SMITH_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 smith::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
mfem::Array< int > dof_list(const fes_t *fes) const
get mfem degree of freedom list for a given FiniteElementSpace
ThermalMaterialIntegrand(MaterialType material)
Construct a ThermalMaterialIntegrand functor with material model of type MaterialType.
auto operator()(double, X x, T temperature, dT_dt dtemp_dt, Params... params) const
Evaluate integrand.
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.