Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
heat_transfer.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
15 #include "mfem.hpp"
16 
17 #include "smith/physics/common.hpp"
20 #include "smith/numerics/odes.hpp"
24 #include "smith/physics/mesh.hpp"
26 
27 namespace smith {
28 
29 namespace heat_transfer {
30 
35  .preconditioner = Preconditioner::HypreL1Jacobi,
36  .relative_tol = 1.0e-6,
37  .absolute_tol = 1.0e-12,
38  .max_iterations = 200};
39 
41 #ifdef MFEM_USE_STRUMPACK
43 #else
45 #endif
46 
51  .relative_tol = 1.0e-4,
52  .absolute_tol = 1.0e-8,
53  .max_iterations = 500,
54  .print_level = 1};
55 
61 
66 
67 } // namespace heat_transfer
68 
81 template <int order, int dim, typename parameters = Parameters<>,
82  typename parameter_indices = std::make_integer_sequence<int, parameters::n>>
84 
86 template <int order, int dim, typename... parameter_space, int... parameter_indices>
87 class HeatTransfer<order, dim, Parameters<parameter_space...>, std::integer_sequence<int, parameter_indices...>>
88  : public BasePhysics {
89  public:
91  static constexpr int VALUE = 0, DERIVATIVE = 1;
92  static constexpr int SHAPE = 0;
93  static constexpr auto I = Identity<dim>();
95 
98  static constexpr auto NUM_STATE_VARS = 2;
99 
117  HeatTransfer(const NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions lin_opts,
118  const smith::TimesteppingOptions timestepping_opts, const std::string& physics_name,
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)
123  {
124  }
125 
142  HeatTransfer(std::unique_ptr<smith::EquationSolver> solver, const smith::TimesteppingOptions timestepping_opts,
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),
147  temperature_(StateManager::newState(H1<order>{}, detail::addPrefix(physics_name, "temperature"), mesh_->tag())),
148  temperature_rate_(
149  StateManager::newState(H1<order>{}, detail::addPrefix(physics_name, "temperature_rate"), mesh_->tag())),
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_)
162  {
163  SLIC_ERROR_ROOT_IF(
164  mfemParMesh().Dimension() != dim,
165  std::format("Compile time class dimension template parameter and runtime mesh dimension do not match"));
166 
167  SLIC_ERROR_ROOT_IF(
168  !nonlin_solver_,
169  "EquationSolver argument is nullptr in HeatTransfer constructor. It is possible that it was previously moved.");
170 
171  // Check for dynamic mode
172  if (timestepping_opts.timestepper != TimestepMethod::QuasiStatic) {
173  ode_.SetTimestepper(timestepping_opts.timestepper);
174  ode_.SetEnforcementMethod(timestepping_opts.enforcement_method);
175  is_quasistatic_ = false;
176  } else {
177  is_quasistatic_ = true;
178  }
179 
180  states_.push_back(&temperature_);
181  if (!is_quasistatic_) {
182  states_.push_back(&temperature_rate_);
183  }
184 
185  adjoints_.push_back(&adjoint_temperature_);
186 
187  // Create a pack of the primal field and parameter finite element spaces
188  const mfem::ParFiniteElementSpace* test_space = &temperature_.space();
189  const mfem::ParFiniteElementSpace* shape_space = &mesh_->shapeDisplacementSpace();
190 
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();
194 
195  SLIC_ERROR_ROOT_IF(
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()));
199 
200  if constexpr (sizeof...(parameter_space) > 0) {
201  tuple<parameter_space...> types{};
202  for_constexpr<sizeof...(parameter_space)>([&](auto i) {
203  parameters_.emplace_back(mfemParMesh(), get<i>(types), detail::addPrefix(name_, parameter_names[i]));
204 
205  trial_spaces[i + NUM_STATE_VARS] = &(parameters_[i].state->space());
206  });
207  }
208 
209  residual_ =
210  std::make_unique<ShapeAwareFunctional<shape_trial, test(scalar_trial, scalar_trial, parameter_space...)>>(
211  shape_space, test_space, trial_spaces);
212 
213  nonlin_solver_->setOperator(residual_with_bcs_);
214 
215  int true_size = temperature_.space().TrueVSize();
216  u_.SetSize(true_size);
217  u_predicted_.SetSize(true_size);
218 
219  initializeThermalStates();
220  }
221 
231  HeatTransfer(const HeatTransferInputOptions& input_options, const std::string& physics_name,
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)
235  {
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));
241  }
242  }
243 
244  if (input_options.initial_temperature) {
245  auto temp = input_options.initial_temperature->constructScalar();
246  temperature_.project(*temp);
247  }
248 
249  if (input_options.source_coef) {
250  // TODO: Not implemented yet in input files
251  // NOTE: cannot use std::functions that use mfem::vector
252  SLIC_ERROR("'source' is not implemented yet in input files.");
253  }
254 
255  // Process the BCs in sorted order for correct behavior with repeated attributes
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) {
259  // FIXME: Better naming for boundary conditions?
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) {
264  // TODO: Not implemented yet in input files
265  // NOTE: cannot use std::functions that use mfem::vector
266  SLIC_ERROR("'flux' is not implemented yet in input files.");
267  } else {
268  SLIC_WARNING_ROOT("Ignoring boundary condition with unknown name: " << physics_name);
269  }
270  }
271  }
272 
277  {
278  dt_ = 0.0;
279  previous_dt_ = -1.0;
280  time_end_step_ = 0.0;
281 
282  u_ = 0.0;
283  temperature_ = 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;
289  }
290 
297  void resetStates(int cycle = 0, double time = 0.0) override
298  {
300  initializeThermalStates();
301 
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));
307  }
308  }
309  }
310 
319  template <typename AppliedTemperatureFunction>
320  void setTemperatureBCs(AppliedTemperatureFunction applied_temperature, Domain& domain)
321  {
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);
325  });
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);
329  }
330 
339  void setTemperatureBCs(const std::set<int>& temp_bdr, std::function<double(const mfem::Vector& x, double t)> temp)
340  {
341  // Project the coefficient onto the grid function
342  temp_bdr_coef_ = std::make_shared<mfem::FunctionCoefficient>(temp);
343 
344  bcs_.addEssential(temp_bdr, temp_bdr_coef_, temperature_.space());
345  }
346 
354  void advanceTimestep(double dt) override
355  {
356  if (is_quasistatic_) {
357  time_ += dt;
358 
359  // Project the essential boundary coefficients
360  for (auto& bc : bcs_.essentials()) {
361  bc.setDofs(temperature_, time_);
362  }
363  nonlin_solver_->solve(temperature_);
364  } else {
365  // Step the time integrator
366  // Note that the ODE solver handles the essential boundary condition application itself
367 
368  // The current ode interface tracks 2 times, one internally which we have a handle to via time_,
369  // and one here via the step interface.
370  // We are ignoring this one, and just using the internal version for now.
371  // This may need to be revisited when more complex time integrators are required,
372  // but at the moment, the double times creates a lot of confusion, so
373  // we short circuit the extra time here by passing a dummy time and ignoring it.
374  double time_tmp = time_;
375  ode_.Step(temperature_, time_tmp, dt);
376  }
377 
378  cycle_ += 1;
379 
380  if (checkpoint_to_disk_) {
381  outputStateToDisk();
382  } else {
383  auto state_names = stateNames();
384  for (const auto& state_name : state_names) {
385  checkpoint_states_[state_name].push_back(state(state_name));
386  }
387  }
388 
389  if (cycle_ > max_cycle_) {
390  timesteps_.push_back(dt);
391  max_cycle_ = cycle_;
392  max_time_ = time_;
393  }
394  }
395 
400  template <typename MaterialType>
401  struct ThermalMaterialIntegrand {
407  ThermalMaterialIntegrand(MaterialType material) : material_(material) {}
408 
412  template <typename X, typename T, typename dT_dt, typename... Params>
413  auto operator()(double /*time*/, X x, T temperature, dT_dt dtemp_dt, Params... params) const
414  {
415  // Get the value and the gradient from the input tuple
416  auto [u, du_dX] = temperature;
417  auto du_dt = get<VALUE>(dtemp_dt);
418 
419  auto [heat_capacity, heat_flux] = material_(x, u, du_dX, params...);
420 
421  return smith::tuple{heat_capacity * du_dt, -1.0 * heat_flux};
422  }
423 
424  private:
425  MaterialType material_;
426  };
427 
452  template <int... active_parameters, typename MaterialType>
453  void setMaterial(DependsOn<active_parameters...>, const MaterialType& material, Domain& domain)
454  {
455  residual_->AddDomainIntegral(Dimension<dim>{}, DependsOn<0, 1, NUM_STATE_VARS + active_parameters...>{},
456  ThermalMaterialIntegrand<MaterialType>(material), domain);
457  }
458 
460  template <typename MaterialType>
461  void setMaterial(const MaterialType& material, Domain& domain)
462  {
463  setMaterial(DependsOn<>{}, material, domain);
464  }
465 
473  void setTemperature(std::function<double(const mfem::Vector& x, double t)> temp)
474  {
475  // Project the coefficient onto the grid function
476  mfem::FunctionCoefficient temp_coef(temp);
477 
478  temp_coef.SetTime(time_);
479  temperature_.project(temp_coef);
480  }
481 
483  void setTemperature(const FiniteElementState temp) { temperature_ = temp; }
484 
507  template <int... active_parameters, typename SourceType>
508  void setSource(DependsOn<active_parameters...>, SourceType source_function, Domain& domain)
509  {
510  residual_->AddDomainIntegral(
511  Dimension<dim>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
512  [source_function](double t, auto x, auto temperature, auto /* dtemp_dt */, auto... params) {
513  // Get the value and the gradient from the input tuple
514  auto [u, du_dX] = temperature;
515 
516  auto source = source_function(x, t, u, du_dX, params...);
517 
518  // Return the source and the flux as a tuple
519  return smith::tuple{-1.0 * source, smith::zero{}};
520  },
521  domain);
522  }
523 
525  template <typename SourceType>
526  void setSource(SourceType source_function, Domain& domain)
527  {
528  setSource(DependsOn<>{}, source_function, domain);
529  }
530 
553  template <int... active_parameters, typename FluxType>
554  void setFluxBCs(DependsOn<active_parameters...>, FluxType flux_function, Domain& domain)
555  {
556  residual_->AddBoundaryIntegral(
557  Dimension<dim - 1>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
558  [flux_function](double t, auto X, auto u, auto /* dtemp_dt */, auto... params) {
559  auto temp = get<VALUE>(u);
560  auto n = cross(get<DERIVATIVE>(X));
561 
562  return flux_function(X, normalize(n), t, temp, params...);
563  },
564  domain);
565  }
566 
568  template <typename FluxType>
569  void setFluxBCs(FluxType flux_function, Domain& domain)
570  {
571  setFluxBCs(DependsOn<>{}, flux_function, domain);
572  }
573 
579  const smith::FiniteElementState& temperature() const { return temperature_; };
580 
586  const smith::FiniteElementState& temperatureRate() const { return temperature_rate_; };
587 
594  const FiniteElementState& state(const std::string& state_name) const override
595  {
596  if (state_name == "temperature") {
597  return temperature_;
598  } else if (state_name == "temperature_rate") {
599  return temperature_rate_;
600  }
601 
602  SLIC_ERROR_ROOT(
603  std::format("State '{}' requested from solid mechanics module '{}', but it doesn't exist", state_name, name_));
604  return temperature_;
605  }
606 
616  void setState(const std::string& state_name, const FiniteElementState& state) override
617  {
618  if (state_name == "temperature") {
619  temperature_ = state;
620  if (!checkpoint_to_disk_) {
621  checkpoint_states_["temperature"][static_cast<size_t>(cycle_)] = temperature_;
622  }
623  return;
624  }
625 
626  SLIC_ERROR_ROOT(
627  std::format("setState for state named '{}' requested from heat transfer module '{}', but it doesn't exist",
628  state_name, name_));
629  }
630 
636  virtual std::vector<std::string> stateNames() const override
637  {
638  if (is_quasistatic_) {
639  return std::vector<std::string>{"temperature"};
640  } else {
641  return std::vector<std::string>{"temperature", "temperature_rate"};
642  }
643  }
644 
668  template <int... active_parameters, typename callable>
670  {
671  residual_->AddBoundaryIntegral(Dimension<dim - 1>{}, DependsOn<0, 1, active_parameters + NUM_STATE_VARS...>{},
672  qfunction, domain);
673  }
674 
681  const FiniteElementState& adjoint(const std::string& state_name) const override
682  {
683  if (state_name == "temperature") {
684  return adjoint_temperature_;
685  }
686 
687  SLIC_ERROR_ROOT(std::format("Adjoint '{}' requested from solid mechanics module '{}', but it doesn't exist",
688  state_name, name_));
689  return adjoint_temperature_;
690  }
691 
697  void completeSetup() override
698  {
699  if (is_quasistatic_) {
700  residual_with_bcs_ = mfem_ext::StdFunctionOperator(
701  temperature_.space().TrueVSize(),
702 
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...);
706 
707  // TODO this copy is required as the sundials solvers do not allow move assignments because of their memory
708  // tracking strategy
709  // See https://github.com/mfem/mfem/issues/3531
710  r = res;
711  r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
712  },
713 
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...);
717  J_ = assemble(drdu);
718  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
719  return *J_;
720  });
721  } else {
722  residual_with_bcs_ = mfem_ext::StdFunctionOperator(
723  temperature_.space().TrueVSize(),
724 
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...);
729 
730  // TODO this copy is required as the sundials solvers do not allow move assignments because of their memory
731  // tracking strategy
732  // See https://github.com/mfem/mfem/issues/3531
733  r = res;
734  r.SetSubVector(bcs_.allEssentialTrueDofs(), 0.0);
735  },
736 
737  [this](const mfem::Vector& du_dt) -> mfem::Operator& {
738  add(1.0, u_, dt_, du_dt, u_predicted_);
739 
740  // K := dR/du
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));
744 
745  // M := dR/du_dot
746  auto M =
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));
750 
751  // J := M + dt K
752  J_.reset(mfem::Add(1.0, *m_mat, dt_, *k_mat));
753  J_e_ = bcs_.eliminateAllEssentialDofsFromMatrix(*J_);
754 
755  return *J_;
756  });
757  }
758 
759  if (checkpoint_to_disk_) {
760  outputStateToDisk();
761  } else {
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));
766  }
767  }
768  }
769 
784  virtual void setAdjointLoad(std::unordered_map<std::string, const smith::FiniteElementDual&> loads) override
785  {
786  SLIC_ERROR_ROOT_IF(loads.size() == 0,
787  "Adjoint load container size must be greater than 0 in the heat transfer module.");
788 
789  auto temp_adjoint_load = loads.find("temperature");
790  auto temp_rate_adjoint_load = loads.find("temperature_rate"); // does not need to be specified
791 
792  SLIC_ERROR_ROOT_IF(temp_adjoint_load == loads.end(), "Adjoint load for \"temperature\" not found.");
793 
794  temperature_adjoint_load_ = temp_adjoint_load->second;
795  // Add the sign correction to move the term to the RHS
796  temperature_adjoint_load_ *= -1.0;
797 
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;
801  }
802  }
803 
834  void reverseAdjointTimestep() override
835  {
836  auto& lin_solver = nonlin_solver_->linearSolver();
837 
838  mfem::HypreParVector adjoint_essential(temperature_adjoint_load_);
839  adjoint_essential = 0.0;
840 
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");
844 
845  cycle_--; // cycle is now at n \in [0,N-1]
846 
847  double dt = getCheckpointedTimestep(cycle_);
848  auto end_step_solution = getCheckpointedStates(cycle_ + 1);
849 
850  temperature_ = end_step_solution.at("temperature");
851 
852  if (is_quasistatic_) {
853  // We store the previous timestep's temperature as the current temperature for use in the lambdas computing the
854  // sensitivities.
855 
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());
860 
861  for (const auto& bc : bcs_.essentials()) {
862  bc.apply(*J_T, temperature_adjoint_load_, adjoint_essential);
863  }
864 
865  lin_solver.SetOperator(*J_T);
866  lin_solver.Mult(temperature_adjoint_load_, adjoint_temperature_);
867  } else {
868  SLIC_ERROR_ROOT_IF(ode_.GetTimestepper() != TimestepMethod::BackwardEuler,
869  "Only backward Euler implemented for transient adjoint heat conduction.");
870 
871  temperature_rate_ = end_step_solution.at("temperature_rate");
872 
873  // K := dR/du
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));
877 
878  // M := dR/du_dot
879  auto M = smith::get<DERIVATIVE>((*residual_)(time_, shapeDisplacement(), temperature_,
880  differentiate_wrt(temperature_rate_),
881  *parameters_[parameter_indices].state...));
882  std::unique_ptr<mfem::HypreParMatrix> m_mat(assemble(M));
883 
884  J_.reset(mfem::Add(1.0, *m_mat, dt, *k_mat));
885  auto J_T = std::unique_ptr<mfem::HypreParMatrix>(J_->Transpose());
886 
887  // recall that temperature_adjoint_load_vector and d_temperature_dt_adjoint_load_vector were already multiplied by
888  // -1 above
889  mfem::HypreParVector modified_RHS(temperature_adjoint_load_);
890  modified_RHS *= dt;
891  modified_RHS.Add(1.0, temperature_rate_adjoint_load_);
892  modified_RHS.Add(-dt, implicit_sensitivity_temperature_start_of_step_);
893 
894  for (const auto& bc : bcs_.essentials()) {
895  bc.apply(*J_T, modified_RHS, adjoint_essential);
896  }
897 
898  lin_solver.SetOperator(*J_T);
899  lin_solver.Mult(modified_RHS, adjoint_temperature_);
900 
901  // This multiply is technically on M transposed. However, this matrix should be symmetric unless
902  // the thermal capacity is a function of the temperature rate of change, which is thermodynamically
903  // impossible, and fortunately not possible with our material interface.
904  // Not doing the transpose here to avoid doing unnecessary work.
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_); // already multiplied by -1
909  }
910 
911  time_end_step_ = time_;
912  time_ -= dt;
913  }
914 
924  FiniteElementDual computeTimestepSensitivity(size_t parameter_field) override
925  {
926  // TODO: the time is likely not being handled correctly on the reverse pass, but we don't
927  // have tests to confirm.
928  auto drdparam = smith::get<DERIVATIVE>(d_residual_d_[parameter_field](time_end_step_));
929  auto drdparam_mat = assemble(drdparam);
930 
931  drdparam_mat->MultTranspose(adjoint_temperature_, *parameters_[parameter_field].sensitivity);
932 
933  return *parameters_[parameter_field].sensitivity;
934  }
935 
945  {
946  auto drdshape =
947  smith::get<DERIVATIVE>((*residual_)(time_end_step_, differentiate_wrt(shapeDisplacement()), temperature_,
948  temperature_rate_, *parameters_[parameter_indices].state...));
949 
950  auto drdshape_mat = assemble(drdshape);
951 
952  drdshape_mat->MultTranspose(adjoint_temperature_, shape_displacement_dual_);
953 
954  return shapeDisplacementSensitivity();
955  }
956 
964  const std::unordered_map<std::string, const smith::FiniteElementDual&> computeInitialConditionSensitivity()
965  const override
966  {
967  return {{"temperature", implicit_sensitivity_temperature_start_of_step_}};
968  }
969 
971  virtual ~HeatTransfer() = default;
972 
973  protected:
976 
979 
981  using test = H1<order>;
982 
985 
988 
991 
994 
997 
1000 
1002  std::unique_ptr<ShapeAwareFunctional<shape_trial, test(scalar_trial, scalar_trial, parameter_space...)>> residual_;
1003 
1005  std::unique_ptr<mfem::HypreParMatrix> M_;
1006 
1008  std::shared_ptr<mfem::Coefficient> temp_bdr_coef_;
1009 
1015 
1017  std::unique_ptr<EquationSolver> nonlin_solver_;
1018 
1025 
1027  std::unique_ptr<mfem::HypreParMatrix> J_;
1028 
1031  std::unique_ptr<mfem::HypreParMatrix> J_e_;
1032 
1034  double dt_;
1035 
1038 
1042 
1044  mfem::Vector u_;
1045 
1047  mfem::Vector u_predicted_;
1048 
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) {
1056  return (*residual_)(DifferentiateWRT<NUM_STATE_VARS + 1 + parameter_indices>{}, TIME, shape_displacement_,
1057  temperature_, temperature_rate_, *parameters_[parameter_indices].state...);
1058  }...};
1059 };
1060 
1061 } // namespace smith
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...
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_ext::FirstOrderODE ode_
the ordinary differential equation that describes how to solve for the time derivative of temperature...
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::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.
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.
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.
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 std::vector< std::string > stateNames() const override
Get a vector of the finite element state primal solution names.
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 ...
Definition: odes.hpp:240
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.
An object containing all input file options for the solver for a heat transfer PDE.
Smith mesh class which assists in constructing the appropriate parallel mfem meshes and registering a...
constexpr SMITH_HOST_DEVICE void for_constexpr(const lambda &f)
multidimensional loop tool that evaluates the lambda body inside the innermost loop.
constexpr auto get(std::integer_sequence< int, n... >)
return the Ith integer in {n...}
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.
Definition: smith.cpp:36
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)
Definition: tensor.hpp:966
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,...
Definition: tensor.hpp:1122
bool holds_alternative(const variant< T0, T1 > &v)
Checks whether a variant's active member is of a certain type.
Definition: variant.hpp:381
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.
Definition: geometry.hpp:17
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
mfem::Array< int > dof_list(const fes_t *fes) const
get mfem degree of freedom list for a given FiniteElementSpace
Definition: domain.cpp:466
H1 elements of order p.
Stores all information held in the input file that is used to configure the solver.
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 ...
Definition: common.hpp:45
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 ...
Definition: tuple.hpp:28
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122
The material and load types for the thermal functional physics module.