Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
odes.cpp
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 
8 
9 #include <utility>
10 #include <vector>
11 
12 namespace smith::mfem_ext {
13 
15  : mfem::SecondOrderTimeDependentOperator(n, 0.0), state_(std::move(state)), solver_(solver), bcs_(bcs), zero_(n)
16 {
17  zero_ = 0.0;
18  U_minus_.SetSize(n);
19  U_.SetSize(n);
20  U_plus_.SetSize(n);
21  dU_dt_.SetSize(n);
22  d2U_dt2_.SetSize(n);
23 }
24 
26 {
27  timestepper_ = timestepper;
28 
29  switch (timestepper) {
31  second_order_ode_solver_ = std::make_unique<mfem::NewmarkSolver>();
32  break;
34  second_order_ode_solver_ = std::make_unique<mfem::HHTAlphaSolver>();
35  break;
37  second_order_ode_solver_ = std::make_unique<mfem::WBZAlphaSolver>();
38  break;
40  // WARNING: apparently mfem's implementation of AverageAccelerationSolver
41  // is NOT equivalent to Newmark (beta = 0.25, gamma = 0.5), so the
42  // stability analysis of using DirichletEnforcementMethod::DirectControl
43  // with Newmark methods does not apply.
44  //
45  // TODO: do a more thorough stability analysis for mfem::GeneralizedAlpha2Solver
46  // to characterize which parameter combinations work with time-varying
47  // dirichlet constraints
48  SLIC_WARNING_ROOT(
49  "Cannot guarantee stability for AverageAcceleration with time-dependent Dirichlet Boundary Conditions");
50  second_order_ode_solver_ = std::make_unique<mfem::AverageAccelerationSolver>();
51  break;
53  second_order_ode_solver_ = std::make_unique<mfem::LinearAccelerationSolver>();
54  break;
56  second_order_ode_solver_ = std::make_unique<mfem::CentralDifferenceSolver>();
57  break;
59  second_order_ode_solver_ = std::make_unique<mfem::FoxGoodwinSolver>();
60  break;
62  first_order_system_ode_solver_ = std::make_unique<mfem::BackwardEulerSolver>();
63  break;
64  default:
65  SLIC_ERROR_ROOT("Timestep method was not a supported second-order ODE method");
66  }
67 
68  if (second_order_ode_solver_) {
69  second_order_ode_solver_->Init(*this);
70  } else if (first_order_system_ode_solver_) {
71  // we need to adjust the width of this operator
72  width *= 2;
73  first_order_system_ode_solver_->Init(*this);
74  } else {
75  // there is no ode_solver
76  SLIC_ERROR("Neither second_order_ode_solver_ nor first_order_system_ode_solver_ specified");
77  }
78 }
79 
80 void SecondOrderODE::Step(mfem::Vector& x, mfem::Vector& dxdt, double& time, double& dt)
81 {
82  if (second_order_ode_solver_) {
83  // if we used a 2nd order method
84  second_order_ode_solver_->Step(x, dxdt, time, dt);
85 
86  if (enforcement_method_ == DirichletEnforcementMethod::FullControl) {
87  U_minus_ = 0.0;
88  U_ = 0.0;
89  U_plus_ = 0.0;
90  for (const auto& bc : bcs_.essentials()) {
91  bc.setDofs(U_minus_, t - epsilon);
92  bc.setDofs(U_, t);
93  bc.setDofs(U_plus_, t + epsilon);
94  }
95 
96  auto constrained_dofs = bcs_.allEssentialTrueDofs();
97  for (int i = 0; i < constrained_dofs.Size(); i++) {
98  x[i] = U_[i];
99  dxdt[i] = (U_plus_[i] - U_minus_[i]) / (2.0 * epsilon);
100  }
101  }
102 
103  } else if (first_order_system_ode_solver_) {
104  // Would be better if displacement and velocity were from a block vector?
105  mfem::Array<int> boffsets(3);
106  boffsets[0] = 0;
107  boffsets[1] = x.Size();
108  boffsets[2] = x.Size() + dxdt.Size();
109  mfem::BlockVector bx(boffsets);
110  bx.GetBlock(0) = x;
111  bx.GetBlock(1) = dxdt;
112 
113  first_order_system_ode_solver_->Step(bx, time, dt);
114 
115  // Copy back
116  x = bx.GetBlock(0);
117  dxdt = bx.GetBlock(1);
118  } else {
119  SLIC_ERROR_ROOT("Neither second_order_ode_solver_ nor first_order_system_ode_solver_ specified");
120  }
121 }
122 
123 void SecondOrderODE::ImplicitSolve(const double dt, const mfem::Vector& u, mfem::Vector& du_dt)
124 {
125  /* A second order o.d.e can be recast as a first order system
126  u_next = u_prev + dt * v_next
127  v_next = v_prev + dt * a_next
128 
129  This means:
130  u_next = u_prev + dt * (v_prev + dt * a_next);
131  u_next = (u_prev + dt * v_prev) + dt*dt*a_next
132  */
133 
134  // Split u in half and du_dt in half
135  mfem::Array<int> boffsets(3);
136  boffsets[0] = 0;
137  boffsets[1] = u.Size() / 2;
138  boffsets[2] = u.Size();
139 
140  const mfem::BlockVector bu(u.GetData(), boffsets);
141 
142  mfem::BlockVector bdu_dt(du_dt.GetData(), boffsets);
143 
144  mfem::Vector u_next(bu.GetBlock(0));
145  u_next.Add(dt, bu.GetBlock(1));
146 
147  Solve(t, dt * dt, dt, u_next,
148  bu.GetBlock(1), // v_next
149  bdu_dt.GetBlock(1)); // a_next
150 
151  bdu_dt.GetBlock(0) = bu.GetBlock(1);
152  bdu_dt.GetBlock(0).Add(dt, bdu_dt.GetBlock(1));
153 }
154 
155 void SecondOrderODE::Solve(const double time, const double c0, const double c1, const mfem::Vector& u,
156  const mfem::Vector& du_dt, mfem::Vector& d2u_dt2) const
157 {
158  // assign these values to variables with greater scope,
159  // so that the residual operator can see them
160  state_.time = time;
161  state_.c0 = c0;
162  state_.c1 = c1;
163  state_.u = u;
164  state_.du_dt = du_dt;
165 
166  // evaluate the constraint functions at a 3-point
167  // stencil of times centered on the time of interest
168  // in order to compute finite-difference approximations
169  // to the time derivatives that appear in the residual
170  U_minus_ = 0.0;
171  U_ = 0.0;
172  U_plus_ = 0.0;
173  for (const auto& bc : bcs_.essentials()) {
174  bc.setDofs(U_minus_, time - epsilon);
175  bc.setDofs(U_, time);
176  bc.setDofs(U_plus_, time + epsilon);
177  }
178 
179  bool implicit = (c0 != 0.0 || c1 != 0.0);
180  if (implicit) {
181  if (enforcement_method_ == DirichletEnforcementMethod::DirectControl) {
182  // TODO: MFEM PR#3064 explicitly deleted the const vector operator-. This
183  // is in active discusion and may be un-deleted. Original line commented.
184  // d2U_dt2_ = (U_ - u) / c0;
185  subtract(1.0 / c0, U_, u, d2U_dt2_);
186  dU_dt_ = du_dt;
187  U_ = u;
188  }
189 
190  if (enforcement_method_ == DirichletEnforcementMethod::RateControl) {
191  // d2U_dt2_ = ((U_plus_ - U_minus_) / (2.0 * epsilon) - du_dt) / c1;
192  subtract(U_plus_, U_minus_, d2U_dt2_);
193  d2U_dt2_ /= 2.0 * epsilon;
194  d2U_dt2_ -= du_dt;
195  d2U_dt2_ /= c1;
196 
197  dU_dt_ = du_dt;
198  U_ = u;
199  }
200 
201  if (enforcement_method_ == DirichletEnforcementMethod::FullControl) {
202  // d2U_dt2_ = (U_minus_ - 2.0 * U_ + U_plus_) / (epsilon * epsilon);
203  add(1.0, U_minus_, -2.0, U_, d2U_dt2_);
204  d2U_dt2_ += U_plus_;
205  d2U_dt2_ /= epsilon * epsilon;
206 
207  // d2U_dt2_ = ((U_plus_ - U_minus_) / (2.0 * epsilon) - du_dt) / c1;
208  subtract(U_plus_, U_minus_, d2U_dt2_);
209  d2U_dt2_ /= 2.0 * epsilon;
210  d2U_dt2_ -= du_dt;
211  d2U_dt2_ /= c1;
212 
213  // dU_dt_ = (U_plus_ - U_minus_) / (2.0 * epsilon) - c1 * d2U_dt2_;
214  subtract(U_plus_, U_minus_, dU_dt_);
215  dU_dt_ /= 2.0 * epsilon;
216  dU_dt_.Add(-1.0 * c1, d2U_dt2_);
217 
218  // U_ = U_ - c0 * d2U_dt2_;
219  U_.Add(-1.0 * c0, d2U_dt2_);
220  }
221  } else {
222  // d2U_dt2_ = (U_minus_ - 2.0 * U_ + U_plus_) / (epsilon * epsilon);
223  add(1.0, U_minus_, -2.0, U_, d2U_dt2_);
224  d2U_dt2_ += U_plus_;
225  d2U_dt2_ /= epsilon * epsilon;
226 
227  // dU_dt_ = (U_plus_ - U_minus_) / (2.0 * epsilon);
228  subtract(U_plus_, U_minus_, dU_dt_);
229  dU_dt_ /= 2.0 * epsilon;
230  }
231 
232  auto constrained_dofs = bcs_.allEssentialTrueDofs();
233  state_.u.SetSubVector(constrained_dofs, 0.0);
234  U_.SetSubVectorComplement(constrained_dofs, 0.0);
235  state_.u += U_;
236 
237  state_.du_dt.SetSubVector(constrained_dofs, 0.0);
238  dU_dt_.SetSubVectorComplement(constrained_dofs, 0.0);
239  state_.du_dt += dU_dt_;
240 
241  // use 0 as our starting guess for unconstrained dofs
242  d2u_dt2 = state_.d2u_dt2;
243  d2u_dt2.SetSubVector(constrained_dofs, 0.0);
244  d2u_dt2 += d2U_dt2_;
245  d2u_dt2.SetSubVectorComplement(constrained_dofs, 0.0);
246 
247  solver_.solve(d2u_dt2);
248  SLIC_WARNING_ROOT_IF(!solver_.nonlinearSolver().GetConverged(), "Newton Solver did not converge.");
249 
250  state_.d2u_dt2 = d2u_dt2;
251 }
252 
254  const BoundaryConditionManager& bcs)
255  : mfem::TimeDependentOperator(n, 0.0), state_(std::move(state)), solver_(solver), bcs_(bcs), zero_(n)
256 {
257  zero_ = 0.0;
258  U_minus_.SetSize(n);
259  U_.SetSize(n);
260  U_plus_.SetSize(n);
261  dU_dt_.SetSize(n);
262 }
263 
265 {
266  timestepper_ = timestepper;
267 
268  switch (timestepper) {
270  ode_solver_ = std::make_unique<mfem::BackwardEulerSolver>();
271  break;
273  ode_solver_ = std::make_unique<mfem::SDIRK33Solver>();
274  break;
276  ode_solver_ = std::make_unique<mfem::ForwardEulerSolver>();
277  break;
279  ode_solver_ = std::make_unique<mfem::RK2Solver>(0.5);
280  break;
282  ode_solver_ = std::make_unique<mfem::RK3SSPSolver>();
283  break;
285  ode_solver_ = std::make_unique<mfem::RK4Solver>();
286  break;
288  ode_solver_ = std::make_unique<mfem::GeneralizedAlphaSolver>(0.5);
289  break;
291  ode_solver_ = std::make_unique<mfem::ImplicitMidpointSolver>();
292  break;
294  ode_solver_ = std::make_unique<mfem::SDIRK23Solver>();
295  break;
297  ode_solver_ = std::make_unique<mfem::SDIRK34Solver>();
298  break;
299  default:
300  SLIC_ERROR_ROOT("Timestep method was not a supported first-order ODE method");
301  }
302  ode_solver_->Init(*this);
303 }
304 
305 void FirstOrderODE::Solve(const double time, const double dt, const mfem::Vector& u, mfem::Vector& du_dt) const
306 {
307  // assign these values to variables with greater scope,
308  // so that the residual operator can see them
309  state_.time = time;
310  state_.dt = dt;
311  state_.u = u;
312 
313  // evaluate the constraint functions at a 3-point
314  // stencil of times centered on the time of interest
315  // in order to compute finite-difference approximations
316  // to the time derivatives that appear in the residual
317  U_minus_ = 0.0;
318  U_ = 0.0;
319  U_plus_ = 0.0;
320  for (const auto& bc : bcs_.essentials()) {
321  bc.setDofs(U_minus_, time - epsilon);
322  bc.setDofs(U_, time);
323  bc.setDofs(U_plus_, time + epsilon);
324  }
325 
326  bool implicit = (dt != 0.0);
327  if (implicit) {
328  if (enforcement_method_ == DirichletEnforcementMethod::DirectControl) {
329  // TODO: MFEM PR#3064 explicitly deleted the const vector operator-. This
330  // is in active discusion and may be un-deleted. Original line commented.
331  // dU_dt_ = (U_ - u) / dt;
332  subtract(1.0 / dt, U_, u, dU_dt_);
333  U_ = u;
334  }
335 
336  if (enforcement_method_ == DirichletEnforcementMethod::RateControl) {
337  // dU_dt_ = (U_plus_ - U_minus_) / (2.0 * epsilon);
338  subtract(1.0 / (2.0 * epsilon), U_plus_, U_minus_, dU_dt_);
339  U_ = u;
340  }
341 
342  if (enforcement_method_ == DirichletEnforcementMethod::FullControl) {
343  // dU_dt_ = (U_plus_ - U_minus_) / (2.0 * epsilon);
344  subtract(1.0 / (2.0 * epsilon), U_plus_, U_minus_, dU_dt_);
345  // U_ = U_ - dt * dU_dt_;
346  U_.Add(-1.0 * dt, dU_dt_);
347  }
348  } else {
349  // dU_dt_ = (U_plus_ - U_minus_) / (2.0 * epsilon);
350  subtract(1.0 / (2.0 * epsilon), U_plus_, U_minus_, dU_dt_);
351  }
352 
353  auto constrained_dofs = bcs_.allEssentialTrueDofs();
354  state_.u.SetSubVector(constrained_dofs, 0.0);
355  U_.SetSubVectorComplement(constrained_dofs, 0.0);
356  state_.u += U_;
357 
358  du_dt = state_.du_dt;
359  du_dt.SetSubVector(constrained_dofs, 0.0);
360  dU_dt_.SetSubVectorComplement(constrained_dofs, 0.0);
361  du_dt += dU_dt_;
362 
363  solver_.solve(du_dt);
364  SLIC_WARNING_ROOT_IF(!solver_.nonlinearSolver().GetConverged(), "Newton Solver did not converge.");
365 
366  state_.du_dt = du_dt;
367  state_.previous_dt = dt;
368 }
369 
370 } // namespace smith::mfem_ext
A container for the boundary condition information relating to a specific physics module.
const mfem::Array< int > & allEssentialTrueDofs() const
Returns all the true degrees of freedom associated with all the essential BCs.
std::vector< BoundaryCondition > & essentials()
Accessor for the essential BC objects.
This class manages the objects typically required to solve a nonlinear set of equations arising from ...
mfem::NewtonSolver & nonlinearSolver()
void solve(mfem::Vector &x) const
FirstOrderODE(int n, FirstOrderODE::State &&state, const EquationSolver &solver, const BoundaryConditionManager &bcs)
Constructor defining the size and specific system of ordinary differential equations to be solved.
Definition: odes.cpp:253
static constexpr double epsilon
a small number used to compute finite difference approximations to time derivatives of boundary condi...
Definition: odes.hpp:251
void SetTimestepper(const smith::TimestepMethod timestepper)
Set the time integration method.
Definition: odes.cpp:264
void ImplicitSolve(const double c0, const double c1, const mfem::Vector &u, const mfem::Vector &du_dt, mfem::Vector &d2u_dt2) override
Solves the equation d2u_dt2 = f(u + c0 * d2u_dt2, du_dt + c1 * d2u_dt2, t)
Definition: odes.hpp:127
SecondOrderODE(int n, State &&state, const EquationSolver &solver, const BoundaryConditionManager &bcs)
Constructor defining the size and specific system of ordinary differential equations to be solved.
Definition: odes.cpp:14
void SetTimestepper(const smith::TimestepMethod timestepper)
Set the time integration method.
Definition: odes.cpp:25
void Step(mfem::Vector &x, mfem::Vector &dxdt, double &time, double &dt)
Performs a time step.
Definition: odes.cpp:80
static constexpr double epsilon
a small number used to compute finite difference approximations to time derivatives of boundary condi...
Definition: odes.hpp:47
TimestepMethod
Timestep method of a solver.
Wrapper classes for using MFEM's ODE solver objects.
A set of references to physics-module-owned variables used by the residual operator.
Definition: odes.hpp:257
double & previous_dt
Previous value of dt.
Definition: odes.hpp:281
mfem::Vector & u
Predicted true DOFs.
Definition: odes.hpp:266
double & time
Time value at which the ODE solver wants to compute a residual.
Definition: odes.hpp:261
mfem::Vector & du_dt
Previous value of du_dt.
Definition: odes.hpp:276
double & dt
Current time step.
Definition: odes.hpp:271
A set of references to physics-module-owned variables used by the residual operator.
Definition: odes.hpp:53
double & c0
coefficient used to calculate updated displacement: u_{n + 1} := u + c0 * d2u_dt2
Definition: odes.hpp:62
mfem::Vector & du_dt
Predicted du_dt.
Definition: odes.hpp:77
mfem::Vector & u
Predicted true DOFs.
Definition: odes.hpp:72
mfem::Vector & d2u_dt2
Previous value of d^2u_dt^2.
Definition: odes.hpp:82
double & c1
coefficient used to calculate updated velocity: du_dt_{n+1} := du_dt + c1 * d2u_dt2
Definition: odes.hpp:67
double & time
Time value at which the ODE solver wants to compute a residual.
Definition: odes.hpp:57