Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_StepperNewmarkImplicitDForm_impl.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_StepperNewmarkImplicitDForm_impl_hpp
11 #define Tempus_StepperNewmarkImplicitDForm_impl_hpp
12 
14 
15 //#define VERBOSE_DEBUG_OUTPUT
16 //#define DEBUG_OUTPUT
17 
18 namespace Tempus {
19 
20 template <class Scalar>
23  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const
24 {
25 #ifdef VERBOSE_DEBUG_OUTPUT
26  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
27 #endif
28  // vPred = v + dt*(1.0-gamma_)*a
29  Thyra::V_StVpStV(Teuchos::ptrFromRef(vPred), 1.0, v, dt * (1.0 - gamma_), a);
30 }
31 
32 template <class Scalar>
36  const Scalar dt) const
37 {
38 #ifdef VERBOSE_DEBUG_OUTPUT
39  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
40 #endif
42  Thyra::createMember<Scalar>(dPred.space());
43  // dPred = dt*v + dt*dt/2.0*(1.0-2.0*beta_)*a
44  Scalar aConst = dt * dt / 2.0 * (1.0 - 2.0 * beta_);
45  Thyra::V_StVpStV(Teuchos::ptrFromRef(dPred), dt, v, aConst, a);
46  // dPred += d;
47  Thyra::Vp_V(Teuchos::ptrFromRef(dPred), d, 1.0);
48 }
49 
50 template <class Scalar>
53  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const
54 {
55 #ifdef VERBOSE_DEBUG_OUTPUT
56  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
57 #endif
58  // v = vPred + dt*gamma_*a
59  Thyra::V_StVpStV(Teuchos::ptrFromRef(v), 1.0, vPred, dt * gamma_, a);
60 }
61 
62 template <class Scalar>
65  const Thyra::VectorBase<Scalar>& a, const Scalar dt) const
66 {
67 #ifdef VERBOSE_DEBUG_OUTPUT
68  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
69 #endif
70  // d = dPred + beta_*dt*dt*a
71  Thyra::V_StVpStV(Teuchos::ptrFromRef(d), 1.0, dPred, beta_ * dt * dt, a);
72 }
73 
74 template <class Scalar>
77  const Thyra::VectorBase<Scalar>& d, const Scalar dt) const
78 {
79 #ifdef VERBOSE_DEBUG_OUTPUT
80  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
81 #endif
82  // a = (d - dPred) / (beta_*dt*dt)
83  Scalar const c = 1.0 / beta_ / dt / dt;
84  Thyra::V_StVpStV(Teuchos::ptrFromRef(a), c, d, -c, dPred);
85 }
86 
87 template <class Scalar>
89 {
90  if (schemeName_ != "User Defined") {
91  *out_ << "\nWARNING: schemeName != 'User Defined' (=" << schemeName_
92  << ").\n"
93  << " Not setting beta, and leaving as beta = " << beta_ << "!\n";
94  return;
95  }
96 
97  beta_ = beta;
98 
99  if (beta_ == 0.0) {
100  *out_ << "\nWARNING: Running (implicit implementation of) Newmark "
101  << "Implicit a-Form Stepper with Beta = 0.0, which \n"
102  << "specifies an explicit scheme. Mass lumping is not possible, "
103  << "so this will be slow! To run explicit \n"
104  << "implementation of Newmark Implicit a-Form Stepper, please "
105  << "re-run with 'Stepper Type' = 'Newmark Explicit a-Form'.\n"
106  << "This stepper allows for mass lumping when called through "
107  << "Piro::TempusSolver.\n";
108  }
109 
111  (beta_ > 1.0) || (beta_ < 0.0), std::logic_error,
112  "\nError in 'Newmark Implicit a-Form' stepper: invalid value of Beta = "
113  << beta_ << ". Please select Beta >= 0 and <= 1. \n");
114 
115  this->isInitialized_ = false;
116 }
117 
118 template <class Scalar>
120 {
121  if (schemeName_ != "User Defined") {
122  *out_ << "\nWARNING: schemeName != 'User Defined' (=" << schemeName_
123  << ").\n"
124  << " Not setting gamma, and leaving as gamma = " << gamma_ << "!\n";
125  return;
126  }
127 
128  gamma_ = gamma;
129 
131  (gamma_ > 1.0) || (gamma_ < 0.0), std::logic_error,
132  "\nError in 'Newmark Implicit a-Form' stepper: invalid value of Gamma ="
133  << gamma_ << ". Please select Gamma >= 0 and <= 1. \n");
134 
135  this->isInitialized_ = false;
136 }
137 
138 template <class Scalar>
140 {
141  schemeName_ = schemeName;
142 
143  if (schemeName_ == "Average Acceleration") {
144  beta_ = 0.25;
145  gamma_ = 0.5;
146  }
147  else if (schemeName_ == "Linear Acceleration") {
148  beta_ = 0.25;
149  gamma_ = 1.0 / 6.0;
150  }
151  else if (schemeName_ == "Central Difference") {
152  beta_ = 0.0;
153  gamma_ = 0.5;
154  }
155  else if (schemeName_ == "User Defined") {
156  beta_ = 0.25;
157  gamma_ = 0.5; // Use defaults until setBeta and setGamma calls.
158  }
159  else {
161  true, std::logic_error,
162  "\nError in Tempus::StepperNewmarkImplicitDForm! "
163  << "Invalid Scheme Name = " << schemeName_ << ". \n"
164  << "Valid Scheme Names are: 'Average Acceleration', "
165  << "'Linear Acceleration', \n"
166  << "'Central Difference' and 'User Defined'.\n");
167  }
168 
169  this->isInitialized_ = false;
170 }
171 
172 template <class Scalar>
174  : out_(Teuchos::VerboseObjectBase::getDefaultOStream())
175 {
176 #ifdef VERBOSE_DEBUG_OUTPUT
177  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
178 #endif
179 
180  this->setStepperName("Newmark Implicit d-Form");
181  this->setStepperType("Newmark Implicit d-Form");
182  this->setUseFSAL(false);
183  this->setICConsistency("None");
184  this->setICConsistencyCheck(false);
185  this->setZeroInitialGuess(false);
186  this->setSchemeName("Average Acceleration");
187 
188  this->setAppAction(Teuchos::null);
189  this->setDefaultSolver();
190 }
191 
192 template <class Scalar>
194  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel,
196  bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck,
197  bool zeroInitialGuess, std::string schemeName, Scalar beta, Scalar gamma,
199  stepperAppAction)
200  : out_(Teuchos::VerboseObjectBase::getDefaultOStream())
201 {
202  this->setStepperName("Newmark Implicit d-Form");
203  this->setStepperType("Newmark Implicit d-Form");
204  this->setUseFSAL(useFSAL);
205  this->setICConsistency(ICConsistency);
206  this->setICConsistencyCheck(ICConsistencyCheck);
207  this->setZeroInitialGuess(zeroInitialGuess);
208  this->setSchemeName(schemeName);
209  this->setBeta(beta);
210  this->setGamma(gamma);
211  this->setAppAction(stepperAppAction);
212  this->setSolver(solver);
213 
214  if (appModel != Teuchos::null) {
215  this->setModel(appModel);
216  this->initialize();
217  }
218 }
219 
220 template <class Scalar>
222  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& appModel)
223 {
224 #ifdef VERBOSE_DEBUG_OUTPUT
225  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
226 #endif
227  validSecondOrderODE_DAE(appModel);
228  this->wrapperModel_ =
230  appModel, "Newmark Implicit d-Form"));
231 
232  TEUCHOS_TEST_FOR_EXCEPTION(this->getSolver() == Teuchos::null,
233  std::logic_error, "Error - Solver is not set!\n");
234  this->getSolver()->setModel(this->wrapperModel_);
235 
236  this->isInitialized_ = false;
237 }
238 
239 template <class Scalar>
242 {
243  if (appAction == Teuchos::null) {
244  // Create default appAction
245  stepperNewmarkImpAppAction_ =
247  }
248  else {
249  stepperNewmarkImpAppAction_ = appAction;
250  }
251 
252  this->isInitialized_ = false;
253 }
254 
255 template <class Scalar>
257  const Teuchos::RCP<SolutionHistory<Scalar>>& solutionHistory)
258 {
259 #ifdef VERBOSE_DEBUG_OUTPUT
260  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
261 #endif
262  this->checkInitialized();
263 
264  using Teuchos::RCP;
265 
266  TEMPUS_FUNC_TIME_MONITOR("Tempus::StepperNewmarkImplicitDForm::takeStep()");
267  {
269  solutionHistory->getNumStates() < 2, std::logic_error,
270  "Error - StepperNewmarkImplicitDForm<Scalar>::takeStep(...)\n"
271  << "Need at least two SolutionStates for NewmarkImplicitDForm.\n"
272  << " Number of States = " << solutionHistory->getNumStates()
273  << "\nTry setting in \"Solution History\" \"Storage Type\" = "
274  << "\"Undo\"\n"
275  << " or \"Storage Type\" = \"Static\" and \"Storage Limit\" = "
276  << "\"2\"\n");
277 
278  auto thisStepper = Teuchos::rcpFromRef(*this);
279  stepperNewmarkImpAppAction_->execute(
280  solutionHistory, thisStepper,
282  Scalar>::ACTION_LOCATION::BEGIN_STEP);
283 
284  RCP<SolutionState<Scalar>> workingState =
285  solutionHistory->getWorkingState();
286  RCP<SolutionState<Scalar>> currentState =
287  solutionHistory->getCurrentState();
288 
290  Teuchos::rcp_dynamic_cast<WrapperModelEvaluatorSecondOrder<Scalar>>(
291  this->wrapperModel_);
292 
293  // Get values of d, v and a from previous step
294  RCP<const Thyra::VectorBase<Scalar>> d_old = currentState->getX();
295  RCP<Thyra::VectorBase<Scalar>> v_old = currentState->getXDot();
296  RCP<Thyra::VectorBase<Scalar>> a_old = currentState->getXDotDot();
297 
298  // Get new values of d, v and a from current workingState
299  //(to be updated here)
300  RCP<Thyra::VectorBase<Scalar>> d_new = workingState->getX();
301  RCP<Thyra::VectorBase<Scalar>> v_new = workingState->getXDot();
302  RCP<Thyra::VectorBase<Scalar>> a_new = workingState->getXDotDot();
303 
304  // Get time and dt
305  const Scalar time = currentState->getTime();
306  const Scalar dt = workingState->getTimeStep();
307  // Update time
308  Scalar t = time + dt;
309 
310 #ifdef DEBUG_OUTPUT
311  Teuchos::Range1D range;
312 
313  *out_ << "\n*** d_old ***\n";
314  RTOpPack::ConstSubVectorView<Scalar> dov;
315  d_old->acquireDetachedView(range, &dov);
316  auto doa = dov.values();
317  for (auto i = 0; i < doa.size(); ++i) *out_ << doa[i] << " ";
318  *out_ << "\n*** d_old ***\n";
319 
320  *out_ << "\n*** v_old ***\n";
321  RTOpPack::ConstSubVectorView<Scalar> vov;
322  v_old->acquireDetachedView(range, &vov);
323  auto voa = vov.values();
324  for (auto i = 0; i < voa.size(); ++i) *out_ << voa[i] << " ";
325  *out_ << "\n*** v_old ***\n";
326 
327  *out_ << "\n*** a_old ***\n";
328  RTOpPack::ConstSubVectorView<Scalar> aov;
329  a_old->acquireDetachedView(range, &aov);
330  auto aoa = aov.values();
331  for (auto i = 0; i < aoa.size(); ++i) *out_ << aoa[i] << " ";
332  *out_ << "\n*** a_old ***\n";
333 #endif
334 
335  // allocate d and v predictors
336  RCP<Thyra::VectorBase<Scalar>> d_pred = Thyra::createMember(d_old->space());
337  RCP<Thyra::VectorBase<Scalar>> v_pred = Thyra::createMember(v_old->space());
338 
339  // compute displacement and velocity predictors
340  predictDisplacement(*d_pred, *d_old, *v_old, *a_old, dt);
341  predictVelocity(*v_pred, *v_old, *a_old, dt);
342 
343 #ifdef DEBUG_OUTPUT
344  *out_ << "\n*** d_pred ***\n";
345  RTOpPack::ConstSubVectorView<Scalar> dpv;
346  d_pred->acquireDetachedView(range, &dpv);
347  auto dpa = dpv.values();
348  for (auto i = 0; i < dpa.size(); ++i) *out_ << dpa[i] << " ";
349  *out_ << "\n*** d_pred ***\n";
350 
351  *out_ << "\n*** v_pred ***\n";
352  RTOpPack::ConstSubVectorView<Scalar> vpv;
353  v_pred->acquireDetachedView(range, &vpv);
354  auto vpa = vpv.values();
355  for (auto i = 0; i < vpa.size(); ++i) *out_ << vpa[i] << " ";
356  *out_ << "\n*** v_pred ***\n";
357 
358 #endif
359  // inject d_pred, v_pred, a and other relevant data into wrapperModel
360  wrapperModel->initializeNewmark(v_pred, d_pred, dt, t, beta_, gamma_);
361 
362  // create initial guess in NOX solver
363  RCP<Thyra::VectorBase<Scalar>> initial_guess =
364  Thyra::createMember(d_pred->space());
365  if ((time == solutionHistory->minTime()) &&
366  (this->initialGuess_ != Teuchos::null)) {
367  // if first time step and initialGuess_ is provided, set initial_guess =
368  // initialGuess_ Throw an exception if initial_guess is not compatible
369  // with solution
370  bool is_compatible =
371  (initial_guess->space())->isCompatible(*this->initialGuess_->space());
373  is_compatible != true, std::logic_error,
374  "Error in Tempus::NemwarkImplicitDForm takeStep(): user-provided "
375  "initial guess'!\n"
376  << "for Newton is not compatible with solution vector!\n");
377  Thyra::copy(*this->initialGuess_, initial_guess.ptr());
378  }
379  else {
380  // Otherwise, set initial guess = diplacement predictor
381  Thyra::copy(*d_pred, initial_guess.ptr());
382  }
383 
384  stepperNewmarkImpAppAction_->execute(
385  solutionHistory, thisStepper,
387  Scalar>::ACTION_LOCATION::BEFORE_SOLVE);
388 
389  // Set d_pred as initial guess for NOX solver, and solve nonlinear system.
390  const Thyra::SolveStatus<Scalar> sStatus =
391  (*(this->solver_)).solve(&*initial_guess);
392 
393  workingState->setSolutionStatus(sStatus); // Converged --> pass.
394 
395  stepperNewmarkImpAppAction_->execute(
396  solutionHistory, thisStepper,
398  Scalar>::ACTION_LOCATION::AFTER_SOLVE);
399 
400  // solve will return converged solution in initial_guess
401  // vector. Copy it here to d_new, to define the new displacement.
402  Thyra::copy(*initial_guess, d_new.ptr());
403 
404  // correct acceleration, velocity
405  correctAcceleration(*a_new, *d_pred, *d_new, dt);
406  correctVelocity(*v_new, *v_pred, *a_new, dt);
407 
408 #ifdef DEBUG_OUTPUT
409  *out_ << "\n*** d_new ***\n";
410  RTOpPack::ConstSubVectorView<Scalar> dnv;
411  d_new->acquireDetachedView(range, &dnv);
412  auto dna = dnv.values();
413  for (auto i = 0; i < dna.size(); ++i) *out_ << dna[i] << " ";
414  *out_ << "\n*** d_new ***\n";
415 
416  *out_ << "\n*** v_new ***\n";
417  RTOpPack::ConstSubVectorView<Scalar> vnv;
418  v_new->acquireDetachedView(range, &vnv);
419  auto vna = vnv.values();
420  for (auto i = 0; i < vna.size(); ++i) *out_ << vna[i] << " ";
421  *out_ << "\n*** v_new ***\n";
422 
423  *out_ << "\n*** a_new ***\n";
424  RTOpPack::ConstSubVectorView<Scalar> anv;
425  a_new->acquireDetachedView(range, &anv);
426  auto ana = anv.values();
427  for (auto i = 0; i < ana.size(); ++i) *out_ << ana[i] << " ";
428  *out_ << "\n*** a_new ***\n";
429 #endif
430 
431  workingState->setOrder(this->getOrder());
432  workingState->computeNorms(currentState);
433 
434  stepperNewmarkImpAppAction_->execute(
435  solutionHistory, thisStepper,
437  Scalar>::ACTION_LOCATION::END_STEP);
438  }
439  return;
440 }
441 
448 template <class Scalar>
451 {
452 #ifdef VERBOSE_DEBUG_OUTPUT
453  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
454 #endif
456  rcp(new StepperState<Scalar>(this->getStepperType()));
457  return stepperState;
458 }
459 
460 template <class Scalar>
462  Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel) const
463 {
464 #ifdef VERBOSE_DEBUG_OUTPUT
465  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
466 #endif
467 
468  out.setOutputToRootOnly(0);
469  out << std::endl;
470  Stepper<Scalar>::describe(out, verbLevel);
471  StepperImplicit<Scalar>::describe(out, verbLevel);
472 
473  out << "--- StepperNewmarkImplicitDForm ---\n";
474  out << " schemeName_ = " << schemeName_ << std::endl;
475  out << " beta_ = " << beta_ << std::endl;
476  out << " gamma_ = " << gamma_ << std::endl;
477  out << "-----------------------------------" << std::endl;
478 }
479 
480 template <class Scalar>
482  Teuchos::FancyOStream& out) const
483 {
484  out.setOutputToRootOnly(0);
485  bool isValidSetup = true;
486 
487  if (!Stepper<Scalar>::isValidSetup(out)) isValidSetup = false;
488 
489  // if ( !StepperImplicit<Scalar>::isValidSetup(out) ) isValidSetup = false;
490  if (this->wrapperModel_->getAppModel() == Teuchos::null) {
491  isValidSetup = false;
492  out << "The application ModelEvaluator is not set!\n";
493  }
494 
495  if (this->wrapperModel_ == Teuchos::null) {
496  isValidSetup = false;
497  out << "The wrapper ModelEvaluator is not set!\n";
498  }
499 
500  if (this->solver_ == Teuchos::null) {
501  isValidSetup = false;
502  out << "The solver is not set!\n";
503  }
504 
505  return isValidSetup;
506 }
507 
508 template <class Scalar>
511 {
512 #ifdef VERBOSE_DEBUG_OUTPUT
513  *out_ << "DEBUG: " << __PRETTY_FUNCTION__ << "\n";
514 #endif
515  auto pl = this->getValidParametersBasicImplicit();
516 
517  auto newmarkPL = Teuchos::parameterList("Newmark Parameters");
518  newmarkPL->set<std::string>("Scheme Name", schemeName_);
519  newmarkPL->set<double>("Beta", beta_);
520  newmarkPL->set<double>("Gamma", gamma_);
521  pl->set("Newmark Parameters", *newmarkPL);
522 
523  return pl;
524 }
525 
526 // Nonmember constructor - ModelEvaluator and ParameterList
527 // ------------------------------------------------------------------------
528 template <class Scalar>
531  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>>& model,
533 {
535  stepper->setStepperImplicitValues(pl);
536 
537  if (pl != Teuchos::null) {
538  if (pl->isSublist("Newmark Parameters")) {
539  auto newmarkPL = pl->sublist("Newmark Parameters", true);
540  std::string schemeName =
541  newmarkPL.get<std::string>("Scheme Name", "Average Acceleration");
542  stepper->setSchemeName(schemeName);
543  if (schemeName == "User Defined") {
544  stepper->setBeta(newmarkPL.get<double>("Beta", 0.25));
545  stepper->setGamma(newmarkPL.get<double>("Gamma", 0.5));
546  }
547  }
548  else {
549  stepper->setSchemeName("Average Acceleration");
550  }
551  }
552 
553  if (model != Teuchos::null) {
554  stepper->setModel(model);
555  stepper->initialize();
556  }
557 
558  return stepper;
559 }
560 
561 } // namespace Tempus
562 #endif // Tempus_StepperNewmarkImplicitDForm_impl_hpp
virtual RCP< const VectorSpaceBase< Scalar > > space() const =0
virtual void setSolver(Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver) override
Set solver.
void correctVelocity(Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
T & get(const std::string &name, T def_value)
Teuchos::RCP< StepperNewmarkImplicitDForm< Scalar > > createStepperNewmarkImplicitDForm(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar >> &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void setStepperName(std::string s)
Set the stepper name.
void correctDisplacement(Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
virtual void initialize()
Initialize after construction and changing input parameters.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
A ModelEvaluator for residual evaluations given a state. This ModelEvaluator takes a state...
Thyra Base interface for time steppers.
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
StepperState is a simple class to hold state information about the stepper.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar >> &solutionHistory)
Take the specified timestep, dt, and return true if successful.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
bool isSublist(const std::string &name) const
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
void predictVelocity(Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const
void setICConsistencyCheck(bool c)
virtual void setAppAction(Teuchos::RCP< StepperNewmarkImplicitDFormAppAction< Scalar >> appAction)
virtual void setZeroInitialGuess(bool zIG)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
Application Action for StepperNewmarkImplicitDForm.
void validSecondOrderODE_DAE(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
virtual void setUseFSAL(bool a)
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar >> &appModel)
void correctAcceleration(Thyra::VectorBase< Scalar > &a, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Scalar dt) const
void setStepperType(std::string s)
Set the stepper type.
void setICConsistency(std::string s)