Tempus  Version of the Day
Time Integration
Tempus_StepperNewmarkImplicitDForm_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_StepperNewmarkImplicitDForm_impl_hpp
10 #define Tempus_StepperNewmarkImplicitDForm_impl_hpp
11 
12 #include "NOX_Thyra.H"
14 #include "Tempus_config.hpp"
15 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
16 
17 //#define VERBOSE_DEBUG_OUTPUT
18 //#define DEBUG_OUTPUT
19 
20 namespace Tempus {
21 
22 // Forward Declaration for recursive includes (this Stepper <--> StepperFactory)
23 template <class Scalar>
24 class StepperFactory;
25 
26 template <class Scalar>
27 void
29  Thyra::VectorBase<Scalar>& vPred, const Thyra::VectorBase<Scalar>& v,
30  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const {
31 #ifdef VERBOSE_DEBUG_OUTPUT
32  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
33 #endif
34  // vPred = v + dt*(1.0-gamma_)*a
35  Thyra::V_StVpStV(Teuchos::ptrFromRef(vPred), 1.0, v, dt * (1.0 - gamma_), a);
36 }
37 
38 template <class Scalar>
39 void
41  Thyra::VectorBase<Scalar>& dPred, const Thyra::VectorBase<Scalar>& d,
42  const Thyra::VectorBase<Scalar>& v, const Thyra::VectorBase<Scalar>& a,
43  const Scalar dt) const {
44 #ifdef VERBOSE_DEBUG_OUTPUT
45  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
46 #endif
47  Teuchos::RCP<const Thyra::VectorBase<Scalar>> tmp =
48  Thyra::createMember<Scalar>(dPred.space());
49  // dPred = dt*v + dt*dt/2.0*(1.0-2.0*beta_)*a
50  Scalar aConst = dt * dt / 2.0 * (1.0 - 2.0 * beta_);
51  Thyra::V_StVpStV(Teuchos::ptrFromRef(dPred), dt, v, aConst, a);
52  // dPred += d;
53  Thyra::Vp_V(Teuchos::ptrFromRef(dPred), d, 1.0);
54 }
55 
56 template <class Scalar>
57 void
59  Thyra::VectorBase<Scalar>& v, const Thyra::VectorBase<Scalar>& vPred,
60  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const {
61 #ifdef VERBOSE_DEBUG_OUTPUT
62  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
63 #endif
64  // v = vPred + dt*gamma_*a
65  Thyra::V_StVpStV(Teuchos::ptrFromRef(v), 1.0, vPred, dt * gamma_, a);
66 }
67 
68 template <class Scalar>
69 void
71  Thyra::VectorBase<Scalar>& d, const Thyra::VectorBase<Scalar>& dPred,
72  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const {
73 #ifdef VERBOSE_DEBUG_OUTPUT
74  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
75 #endif
76  // d = dPred + beta_*dt*dt*a
77  Thyra::V_StVpStV(Teuchos::ptrFromRef(d), 1.0, dPred, beta_ * dt * dt, a);
78 }
79 
80 template <class Scalar>
81 void
83  Thyra::VectorBase<Scalar>& a, const Thyra::VectorBase<Scalar>& dPred,
84  const Thyra::VectorBase<Scalar>& d, const Scalar dt) const {
85 #ifdef VERBOSE_DEBUG_OUTPUT
86  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
87 #endif
88  // a = (d - dPred) / (beta_*dt*dt)
89  Scalar const c = 1.0 / beta_ / dt / dt;
90  Thyra::V_StVpStV(Teuchos::ptrFromRef(a), c, d, -c, dPred);
91 }
92 
93 // StepperNewmarkImplicitDForm definitions:
94 template <class Scalar>
96  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel,
97  Teuchos::RCP<Teuchos::ParameterList> pList)
98  : out_(Teuchos::VerboseObjectBase::getDefaultOStream()) {
99 #ifdef VERBOSE_DEBUG_OUTPUT
100  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
101 #endif
102  using Teuchos::RCP;
103  using Teuchos::ParameterList;
104 
105  // Set all the input parameters and call initialize
106  this->setParameterList(pList);
107  this->setModel(appModel);
108  this->initialize();
109 }
110 
111 template <class Scalar>
112 void
114  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel) {
115 #ifdef VERBOSE_DEBUG_OUTPUT
116  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
117 #endif
118  this->validSecondOrderODE_DAE(appModel);
119  Teuchos::RCP<WrapperModelEvaluatorSecondOrder<Scalar> > wrapperModel =
121  appModel, "Newmark Implicit d-Form"));
122  this->wrapperModel_ = wrapperModel;
123  inArgs_ = this->wrapperModel_->getNominalValues();
124  outArgs_ = this->wrapperModel_->createOutArgs();
125 }
126 
127 template <class Scalar>
128 void
130 {
131  TEUCHOS_TEST_FOR_EXCEPTION( this->wrapperModel_ == Teuchos::null,
132  std::logic_error,
133  "Error - Need to set the model, setModel(), before calling "
134  "StepperNewmarkImplicitDForm::initialize()\n");
135 
136 #ifdef VERBOSE_DEBUG_OUTPUT
137  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
138 #endif
139  this->setParameterList(this->stepperPL_);
140  this->setSolver();
141 }
142 
143 template <class Scalar>
144 void
146  const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory) {
147 #ifdef VERBOSE_DEBUG_OUTPUT
148  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
149 #endif
150  using Teuchos::RCP;
151 
152  TEMPUS_FUNC_TIME_MONITOR("Tempus::StepperNewmarkImplicitDForm::takeStep()");
153  {
154  TEUCHOS_TEST_FOR_EXCEPTION(solutionHistory->getNumStates() < 2,
155  std::logic_error,
156  "Error - StepperNewmarkImplicitDForm<Scalar>::takeStep(...)\n"
157  "Need at least two SolutionStates for NewmarkImplicitDForm.\n"
158  " Number of States = " << solutionHistory->getNumStates() << "\n"
159  "Try setting in \"Solution History\" \"Storage Type\" = \"Undo\"\n"
160  " or \"Storage Type\" = \"Static\" and \"Storage Limit\" = \"2\"\n");
161 
162  RCP<SolutionState<Scalar>> workingState =solutionHistory->getWorkingState();
163  RCP<SolutionState<Scalar>> currentState =solutionHistory->getCurrentState();
164 
165  Teuchos::RCP<WrapperModelEvaluatorSecondOrder<Scalar> > wrapperModel =
166  Teuchos::rcp_dynamic_cast<WrapperModelEvaluatorSecondOrder<Scalar> >(
167  this->wrapperModel_);
168 
169  // Get values of d, v and a from previous step
170  RCP<const Thyra::VectorBase<Scalar>> d_old = currentState->getX();
171  RCP<Thyra::VectorBase<Scalar>> v_old = currentState->getXDot();
172  RCP<Thyra::VectorBase<Scalar>> a_old = currentState->getXDotDot();
173 
174  // Get new values of d, v and a from current workingState
175  //(to be updated here)
176  RCP<Thyra::VectorBase<Scalar>> d_new = workingState->getX();
177  RCP<Thyra::VectorBase<Scalar>> v_new = workingState->getXDot();
178  RCP<Thyra::VectorBase<Scalar>> a_new = workingState->getXDotDot();
179 
180  // Get time and dt
181  const Scalar time = currentState->getTime();
182  const Scalar dt = workingState->getTimeStep();
183  // Update time
184  Scalar t = time + dt;
185 
186 #if 0
187  // Compute initial displacement, velocity and acceleration
188  if (time == solutionHistory->minTime()) {
189 
190  RCP<Thyra::VectorBase<Scalar>>
191  d_init = Thyra::createMember(d_old->space());
192 
193  RCP<Thyra::VectorBase<Scalar>>
194  v_init = Thyra::createMember(v_old->space());
195 
196  RCP<Thyra::VectorBase<Scalar>>
197  a_init = Thyra::createMember(a_old->space());
198 
199  Thyra::copy(*a_old, a_init.ptr());
200  Thyra::copy(*v_old, v_init.ptr());
201  Thyra::copy(*d_old, d_init.ptr());
202 
203 #ifdef DEBUG_OUTPUT
204  Teuchos::Range1D range;
205 
206  *out_ << "\n*** d_init ***\n";
207  RTOpPack::ConstSubVectorView<Scalar> div;
208  d_init->acquireDetachedView(range, &div);
209  auto dia = div.values();
210  for (auto i = 0; i < dia.size(); ++i) *out_ << dia[i] << " ";
211  *out_ << "\n*** d_init ***\n";
212 
213  *out_ << "\n*** v_init ***\n";
214  RTOpPack::ConstSubVectorView<Scalar> viv;
215  v_init->acquireDetachedView(range, &viv);
216  auto via = viv.values();
217  for (auto i = 0; i < via.size(); ++i) *out_ << via[i] << " ";
218  *out_ << "\n*** v_init ***\n";
219 
220  *out_ << "\n*** a_init ***\n";
221  RTOpPack::ConstSubVectorView<Scalar> aiv;
222  a_init->acquireDetachedView(range, &aiv);
223  auto aia = aiv.values();
224  for (auto i = 0; i < aia.size(); ++i) *out_ << aia[i] << " ";
225  *out_ << "\n*** a_init ***\n";
226 #endif
227 
228  wrapperModel->initializeNewmark(
229  a_init, v_init, d_init, dt, time, beta_, gamma_);
230 
231  const Thyra::SolveStatus<Scalar>
232  status = this->solveNonLinear(this->wrapperModel_, *this->solver_, d_init, inArgs_);
233 
234  if (status.solveStatus == Thyra::SOLVE_STATUS_CONVERGED ) {
235  workingState->setSolutionStatus(Status::PASSED);
236  } else {
237  workingState->setSolutionStatus(Status::FAILED);
238  }
239 
240  correctAcceleration(*a_old, *d_old, *d_init, dt);
241  Thyra::copy(*d_init, d_old.ptr());
242  }
243 #endif
244 
245 #ifdef DEBUG_OUTPUT
246  Teuchos::Range1D range;
247 
248  *out_ << "\n*** d_old ***\n";
249  RTOpPack::ConstSubVectorView<Scalar> dov;
250  d_old->acquireDetachedView(range, &dov);
251  auto doa = dov.values();
252  for (auto i = 0; i < doa.size(); ++i) *out_ << doa[i] << " ";
253  *out_ << "\n*** d_old ***\n";
254 
255  *out_ << "\n*** v_old ***\n";
256  RTOpPack::ConstSubVectorView<Scalar> vov;
257  v_old->acquireDetachedView(range, &vov);
258  auto voa = vov.values();
259  for (auto i = 0; i < voa.size(); ++i) *out_ << voa[i] << " ";
260  *out_ << "\n*** v_old ***\n";
261 
262  *out_ << "\n*** a_old ***\n";
263  RTOpPack::ConstSubVectorView<Scalar> aov;
264  a_old->acquireDetachedView(range, &aov);
265  auto aoa = aov.values();
266  for (auto i = 0; i < aoa.size(); ++i) *out_ << aoa[i] << " ";
267  *out_ << "\n*** a_old ***\n";
268 #endif
269 
270  // allocate d and v predictors
271  RCP<Thyra::VectorBase<Scalar>> d_pred = Thyra::createMember(d_old->space());
272  RCP<Thyra::VectorBase<Scalar>> v_pred = Thyra::createMember(v_old->space());
273 
274  // compute displacement and velocity predictors
275  predictDisplacement(*d_pred, *d_old, *v_old, *a_old, dt);
276  predictVelocity(*v_pred, *v_old, *a_old, dt);
277 
278 #ifdef DEBUG_OUTPUT
279  *out_ << "\n*** d_pred ***\n";
280  RTOpPack::ConstSubVectorView<Scalar> dpv;
281  d_pred->acquireDetachedView(range, &dpv);
282  auto dpa = dpv.values();
283  for (auto i = 0; i < dpa.size(); ++i) *out_ << dpa[i] << " ";
284  *out_ << "\n*** d_pred ***\n";
285 
286  *out_ << "\n*** v_pred ***\n";
287  RTOpPack::ConstSubVectorView<Scalar> vpv;
288  v_pred->acquireDetachedView(range, &vpv);
289  auto vpa = vpv.values();
290  for (auto i = 0; i < vpa.size(); ++i) *out_ << vpa[i] << " ";
291  *out_ << "\n*** v_pred ***\n";
292 
293 #endif
294  // inject d_pred, v_pred, a and other relevant data into wrapperModel
295  wrapperModel->initializeNewmark(v_pred, d_pred, dt, t, beta_, gamma_);
296 
297  // create initial guess in NOX solver
298  RCP<Thyra::VectorBase<Scalar>> initial_guess = Thyra::createMember(d_pred->space());
299  if ((time == solutionHistory->minTime()) && (initial_guess_ != Teuchos::null)) {
300  //if first time step and initial_guess_ is provided, set initial_guess = initial_guess_
301  //Throw an exception if initial_guess is not compatible with solution
302  bool is_compatible = (initial_guess->space())->isCompatible(*initial_guess_->space());
303  TEUCHOS_TEST_FOR_EXCEPTION(
304  is_compatible != true, std::logic_error,
305  "Error in Tempus::NemwarkImplicitDForm takeStep(): user-provided initial guess'!\n"
306  << "for Newton is not compatible with solution vector!\n");
307  Thyra::copy(*initial_guess_, initial_guess.ptr());
308  }
309  else {
310  //Otherwise, set initial guess = diplacement predictor
311  Thyra::copy(*d_pred, initial_guess.ptr());
312  }
313 
314  //Set d_pred as initial guess for NOX solver, and solve nonlinear system.
315  const Thyra::SolveStatus<Scalar> status =
316  this->solveImplicitODE(initial_guess);
317 
318  if (status.solveStatus == Thyra::SOLVE_STATUS_CONVERGED)
319  workingState->setSolutionStatus(Status::PASSED);
320  else
321  workingState->setSolutionStatus(Status::FAILED);
322 
323  //solveImplicitODE will return converged solution in initial_guess
324  //vector. Copy it here to d_new, to define the new displacement.
325  Thyra::copy(*initial_guess, d_new.ptr());
326 
327  //correct acceleration, velocity
328  correctAcceleration(*a_new, *d_pred, *d_new, dt);
329  correctVelocity(*v_new, *v_pred, *a_new, dt);
330 
331 #ifdef DEBUG_OUTPUT
332  *out_ << "\n*** d_new ***\n";
333  RTOpPack::ConstSubVectorView<Scalar> dnv;
334  d_new->acquireDetachedView(range, &dnv);
335  auto dna = dnv.values();
336  for (auto i = 0; i < dna.size(); ++i) *out_ << dna[i] << " ";
337  *out_ << "\n*** d_new ***\n";
338 
339  *out_ << "\n*** v_new ***\n";
340  RTOpPack::ConstSubVectorView<Scalar> vnv;
341  v_new->acquireDetachedView(range, &vnv);
342  auto vna = vnv.values();
343  for (auto i = 0; i < vna.size(); ++i) *out_ << vna[i] << " ";
344  *out_ << "\n*** v_new ***\n";
345 
346  *out_ << "\n*** a_new ***\n";
347  RTOpPack::ConstSubVectorView<Scalar> anv;
348  a_new->acquireDetachedView(range, &anv);
349  auto ana = anv.values();
350  for (auto i = 0; i < ana.size(); ++i) *out_ << ana[i] << " ";
351  *out_ << "\n*** a_new ***\n";
352 #endif
353 
354  workingState->setOrder(this->getOrder());
355  }
356  return;
357 }
358 
359 /** \brief Provide a StepperState to the SolutionState.
360  * This Stepper does not have any special state data,
361  * so just provide the base class StepperState with the
362  * Stepper description. This can be checked to ensure
363  * that the input StepperState can be used by this Stepper.
364  */
365 template <class Scalar>
366 Teuchos::RCP<Tempus::StepperState<Scalar>>
368 #ifdef VERBOSE_DEBUG_OUTPUT
369  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
370 #endif
371  Teuchos::RCP<Tempus::StepperState<Scalar>> stepperState =
372  rcp(new StepperState<Scalar>(description()));
373  return stepperState;
374 }
375 
376 template <class Scalar>
377 std::string
379 #ifdef VERBOSE_DEBUG_OUTPUT
380  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
381 #endif
382  std::string name = "Newmark Implicit d-Form";
383  return (name);
384 }
385 
386 template <class Scalar>
387 void
389  Teuchos::FancyOStream& out,
390  const Teuchos::EVerbosityLevel verbLevel) const {
391 #ifdef VERBOSE_DEBUG_OUTPUT
392  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
393 #endif
394  out << description() << "::describe:" << std::endl
395  << "wrapperModel_ = " << this->wrapperModel_->description() << std::endl;
396 }
397 
398 template <class Scalar>
399 void
401  Teuchos::RCP<Teuchos::ParameterList> const& pList) {
402 #ifdef VERBOSE_DEBUG_OUTPUT
403  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
404 #endif
405  if (pList == Teuchos::null) {
406  // Create default parameters if null, otherwise keep current parameters.
407  if (this->stepperPL_ == Teuchos::null) this->stepperPL_ = this->getDefaultParameters();
408  } else {
409  this->stepperPL_ = pList;
410  }
411  // Can not validate because of optional Parameters.
412  // stepperPL_->validateParametersAndSetDefaults(*this->getValidParameters());
413  // Get beta and gamma from parameter list
414  // IKT, FIXME: does parameter list get validated somewhere?
415  // validateParameters above is commented out...
416 
417  Teuchos::RCP<Teuchos::ParameterList> stepperPL = this->stepperPL_;
418  std::string stepperType = stepperPL->get<std::string>("Stepper Type");
419  TEUCHOS_TEST_FOR_EXCEPTION(
420  stepperType != "Newmark Implicit d-Form", std::logic_error,
421  "Error - Stepper Type is not 'Newmark Implicit d-Form'!\n"
422  << " Stepper Type = " << stepperPL->get<std::string>("Stepper Type")
423  << "\n");
424  beta_ = 0.25; // default value
425  gamma_ = 0.5; // default value
426  Teuchos::VerboseObjectBase::getDefaultOStream();
427  if (this->stepperPL_->isSublist("Newmark Parameters")) {
428  Teuchos::ParameterList& newmarkPL =
429  this->stepperPL_->sublist("Newmark Parameters", true);
430  std::string scheme_name = newmarkPL.get("Scheme Name", "Not Specified");
431  if (scheme_name == "Not Specified") {
432  beta_ = newmarkPL.get("Beta", 0.25);
433  gamma_ = newmarkPL.get("Gamma", 0.5);
434  TEUCHOS_TEST_FOR_EXCEPTION(
435  (beta_ > 1.0) || (beta_ < 0.0), std::logic_error,
436  "\nError in 'Newmark Implicit d-Form' stepper: invalid value of Beta = "
437  << beta_ << ". Please select Beta >= 0 and <= 1. \n");
438  TEUCHOS_TEST_FOR_EXCEPTION(
439  (gamma_ > 1.0) || (gamma_ < 0.0), std::logic_error,
440  "\nError in 'Newmark Implicit d-Form' stepper: invalid value of Gamma = "
441  << gamma_ << ". Please select Gamma >= 0 and <= 1. \n");
442  *out_ << "\nSetting Beta = " << beta_ << " and Gamma = " << gamma_
443  << " from Newmark Parameters in input file.\n";
444  } else {
445  *out_ << "\nScheme Name = " << scheme_name << ". Using values \n"
446  << "of Beta and Gamma for this scheme (ignoring values of Beta and "
447  "Gamma \n"
448  << "in input file, if provided).\n";
449  if (scheme_name == "Average Acceleration") {
450  beta_ = 0.25;
451  gamma_ = 0.5;
452  } else if (scheme_name == "Linear Acceleration") {
453  beta_ = 0.25;
454  gamma_ = 1.0 / 6.0;
455  } else if (scheme_name == "Central Difference") {
456  beta_ = 0.0;
457  gamma_ = 0.5;
458  } else {
459  TEUCHOS_TEST_FOR_EXCEPTION(
460  true, std::logic_error,
461  "\nError in Tempus::StepperNewmarkImplicitDForm! Invalid Scheme Name = "
462  << scheme_name << ". \n"
463  << "Valid Scheme Names are: 'Average Acceleration', 'Linear "
464  "Acceleration', \n"
465  << "'Central Difference' and 'Not Specified'.\n");
466  }
467  *out_ << "===> Beta = " << beta_ << ", Gamma = " << gamma_ << "\n";
468  }
469  if (beta_ == 0.0) {
470  //IKT, FIXME - THROW EXCEPTION.
471  TEUCHOS_TEST_FOR_EXCEPTION(
472  true, std::logic_error,
473  "Error - d-Form of Newmark scheme is not defined for explicit (Beta = 0).\n");
474  }
475  } else {
476  *out_ << "\nNo Newmark Parameters sublist found in input file; using "
477  "default values of Beta = "
478  << beta_ << " and Gamma = " << gamma_ << ".\n";
479  }
480 }
481 
482 template <class Scalar>
483 Teuchos::RCP<const Teuchos::ParameterList>
485 #ifdef VERBOSE_DEBUG_OUTPUT
486  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
487 #endif
488  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
489  pl->setName("Default Stepper - " + this->description());
490  pl->set("Stepper Type", this->description());
491  pl->set("Zero Initial Guess", false);
492  pl->set("Solver Name", "",
493  "Name of ParameterList containing the solver specifications.");
494 
495  return pl;
496 }
497 template <class Scalar>
498 Teuchos::RCP<Teuchos::ParameterList>
500 #ifdef VERBOSE_DEBUG_OUTPUT
501  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
502 #endif
503  using Teuchos::RCP;
504  using Teuchos::ParameterList;
505 
506  RCP<ParameterList> pl = Teuchos::parameterList();
507  pl->setName("Default Stepper - " + this->description());
508  pl->set<std::string>("Stepper Type", this->description());
509  pl->set<bool> ("Zero Initial Guess", false);
510  pl->set<std::string>("Solver Name", "Default Solver");
511 
512  RCP<ParameterList> solverPL = this->defaultSolverParameters();
513  pl->set("Default Solver", *solverPL);
514 
515  return pl;
516 }
517 
518 template <class Scalar>
519 Teuchos::RCP<Teuchos::ParameterList>
521 #ifdef VERBOSE_DEBUG_OUTPUT
522  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
523 #endif
524  return (this->stepperPL_);
525 }
526 
527 template <class Scalar>
528 Teuchos::RCP<Teuchos::ParameterList>
530 #ifdef VERBOSE_DEBUG_OUTPUT
531  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
532 #endif
533  Teuchos::RCP<Teuchos::ParameterList> temp_plist = this->stepperPL_;
534  this->stepperPL_ = Teuchos::null;
535  return (temp_plist);
536 }
537 
538 } // namespace Tempus
539 #endif // Tempus_StepperNewmarkImplicitDForm_impl_hpp
Tempus::StepperState
StepperState is a simple class to hold state information about the stepper.
Definition: Tempus_StepperState.hpp:36
Tempus::StepperNewmarkImplicitDForm::StepperNewmarkImplicitDForm
StepperNewmarkImplicitDForm()
Default Constructor – not allowed.
Tempus::solutionHistory
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
Definition: Tempus_SolutionHistory_impl.hpp:504
Tempus::StepperNewmarkImplicitDForm::correctDisplacement
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:70
Tempus
Definition: Tempus_AdjointAuxSensitivityModelEvaluator_decl.hpp:20
Tempus::StepperNewmarkImplicitDForm::setParameterList
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:400
Tempus::StepperNewmarkImplicitDForm::getValidParameters
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:484
Tempus::StepperNewmarkImplicitDForm::description
virtual std::string description() const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:378
Tempus::StepperNewmarkImplicitDForm::correctAcceleration
void correctAcceleration(Thyra::VectorBase< Scalar > &a, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Scalar dt) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:82
Tempus::StepperNewmarkImplicitDForm::predictDisplacement
void predictDisplacement(Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:40
Tempus::StepperNewmarkImplicitDForm::takeStep
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar >> &solutionHistory)
Take the specified timestep, dt, and return true if successful.
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:145
Tempus::StepperNewmarkImplicitDForm::getDefaultParameters
Teuchos::RCP< Teuchos::ParameterList > getDefaultParameters() const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:499
Tempus::WrapperModelEvaluatorSecondOrder
A ModelEvaluator for residual evaluations given a state. This ModelEvaluator takes a state,...
Definition: Tempus_WrapperModelEvaluatorSecondOrder_decl.hpp:32
Tempus::StepperNewmarkImplicitDForm::out_
Teuchos::RCP< Teuchos::FancyOStream > out_
Definition: Tempus_StepperNewmarkImplicitDForm_decl.hpp:154
Tempus::StepperNewmarkImplicitDForm::setModel
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar >> &appModel)
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:113
Tempus::PASSED
Definition: Tempus_Types.hpp:17
Tempus::StepperNewmarkImplicitDForm::describe
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:388
Tempus::SolutionHistory
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Definition: Tempus_Integrator.hpp:25
Tempus::StepperNewmarkImplicitDForm::unsetParameterList
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:529
Tempus::FAILED
Definition: Tempus_Types.hpp:18
Tempus::StepperNewmarkImplicitDForm::initialize
virtual void initialize()
Initialize during construction and after changing input parameters.
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:129
Tempus::StepperNewmarkImplicitDForm::getNonconstParameterList
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:520
Tempus::WrapperModelEvaluatorSecondOrder::initializeNewmark
void initializeNewmark(Teuchos::RCP< Vector > v_pred, Teuchos::RCP< Vector > d_pred, Scalar delta_t, Scalar t, Scalar beta, Scalar gamma)
Set values needed in evalModelImpl.
Definition: Tempus_WrapperModelEvaluatorSecondOrder_decl.hpp:81
Tempus::StepperNewmarkImplicitDForm::correctVelocity
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:58
Tempus_StepperFactory.hpp
Teuchos
Definition: Tempus_Integrator.hpp:19
Tempus::StepperNewmarkImplicitDForm::getDefaultStepperState
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:367
Tempus::StepperNewmarkImplicitDForm::predictVelocity
void predictVelocity(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
Definition: Tempus_StepperNewmarkImplicitDForm_impl.hpp:28