9 #ifndef Tempus_StepperIMEX_RK_Partition_impl_hpp
10 #define Tempus_StepperIMEX_RK_Partition_impl_hpp
12 #include "Tempus_config.hpp"
14 #include "Tempus_WrapperModelEvaluatorPairPartIMEX_Basic.hpp"
15 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
16 #include "Thyra_VectorStdOps.hpp"
17 #include "NOX_Thyra.H"
23 template<
class Scalar>
class StepperFactory;
26 template<
class Scalar>
29 this->setStepperType(
"Partitioned IMEX RK SSP2");
30 this->setUseFSAL( this->getUseFSALDefault());
31 this->setICConsistency( this->getICConsistencyDefault());
32 this->setICConsistencyCheck( this->getICConsistencyCheckDefault());
33 this->setZeroInitialGuess(
false);
35 this->setStageNumber(-1);
37 this->setTableaus(
"Partitioned IMEX RK SSP2");
38 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
41 this->setAppAction(Teuchos::null);
42 this->setDefaultSolver();
46 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
47 template<
class Scalar>
49 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
51 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
53 std::string ICConsistency,
54 bool ICConsistencyCheck,
55 bool zeroInitialGuess,
56 std::string stepperType,
61 this->setStepperType( stepperType);
62 this->setUseFSAL( useFSAL);
63 this->setICConsistency( ICConsistency);
64 this->setICConsistencyCheck( ICConsistencyCheck);
65 this->setZeroInitialGuess( zeroInitialGuess);
67 this->setStageNumber(-1);
69 this->setExplicitTableau(explicitTableau);
70 this->setImplicitTableau(implicitTableau);
71 this->setObserver(obs);
72 this->setAppAction(Teuchos::null);
73 this->setSolver(solver);
75 if (appModel != Teuchos::null) {
76 this->setModel(appModel);
81 template<
class Scalar>
83 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel,
84 const Teuchos::RCP<Thyra::NonlinearSolverBase<Scalar> >& solver,
86 std::string ICConsistency,
87 bool ICConsistencyCheck,
88 bool zeroInitialGuess,
90 std::string stepperType,
95 this->setStepperType( stepperType);
96 this->setUseFSAL( useFSAL);
97 this->setICConsistency( ICConsistency);
98 this->setICConsistencyCheck( ICConsistencyCheck);
99 this->setZeroInitialGuess( zeroInitialGuess);
100 this->setOrder( order);
102 this->setStageNumber(-1);
104 this->setExplicitTableau(explicitTableau);
105 this->setImplicitTableau(implicitTableau);
106 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
107 this->setObserver(Teuchos::null);
109 this->setAppAction(stepperRKAppAction);
110 this->setSolver(solver);
112 if (appModel != Teuchos::null) {
113 this->setModel(appModel);
119 template<
class Scalar>
121 std::string stepperType,
125 if (stepperType ==
"") stepperType =
"Partitioned IMEX RK SSP2";
127 if (stepperType ==
"Partitioned IMEX RK 1st order") {
130 typedef Teuchos::ScalarTraits<Scalar> ST;
132 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
133 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
134 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
135 const Scalar one = ST::one();
136 const Scalar zero = ST::zero();
139 A(0,0) = zero; A(0,1) = zero;
140 A(1,0) = one; A(1,1) = zero;
143 b(0) = one; b(1) = zero;
146 c(0) = zero; c(1) = one;
151 "Explicit Tableau - Partitioned IMEX RK 1st order",
152 A,b,c,order,order,order));
154 this->setExplicitTableau(expTableau);
158 typedef Teuchos::ScalarTraits<Scalar> ST;
160 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
161 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
162 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
163 const Scalar one = ST::one();
164 const Scalar zero = ST::zero();
167 A(0,0) = zero; A(0,1) = zero;
168 A(1,0) = zero; A(1,1) = one;
171 b(0) = zero; b(1) = one;
174 c(0) = zero; c(1) = one;
179 "Implicit Tableau - Partitioned IMEX RK 1st order",
180 A,b,c,order,order,order));
182 this->setImplicitTableau(impTableau);
184 this->setStepperType(
"Partitioned IMEX RK 1st order");
187 }
else if (stepperType ==
"Partitioned IMEX RK SSP2") {
190 this->setExplicitTableau(stepperERK->getTableau());
194 stepperSDIRK->setGammaType(
"2nd Order L-stable");
195 this->setImplicitTableau(stepperSDIRK->getTableau());
197 this->setStepperType(
"Partitioned IMEX RK SSP2");
200 }
else if (stepperType ==
"Partitioned IMEX RK ARS 233") {
201 typedef Teuchos::ScalarTraits<Scalar> ST;
203 Teuchos::SerialDenseMatrix<int,Scalar> A(NumStages,NumStages);
204 Teuchos::SerialDenseVector<int,Scalar> b(NumStages);
205 Teuchos::SerialDenseVector<int,Scalar> c(NumStages);
206 const Scalar one = ST::one();
207 const Scalar zero = ST::zero();
208 const Scalar onehalf = ST::one()/(2*ST::one());
209 const Scalar gamma = (3*one+ST::squareroot(3*one))/(6*one);
213 A(0,0) = zero; A(0,1) = zero; A(0,2) = zero;
214 A(1,0) = gamma; A(1,1) = zero; A(1,2) = zero;
215 A(2,0) = (gamma-1.0); A(2,1) = (2.0-2.0*gamma); A(2,2) = zero;
218 b(0) = zero; b(1) = onehalf; b(2) = onehalf;
221 c(0) = zero; c(1) = gamma; c(2) = one-gamma;
226 "Explicit Tableau - Partitioned IMEX RK ARS 233",
227 A,b,c,order,order,order));
229 this->setExplicitTableau(expTableau);
234 A(0,0) = zero; A(0,1) = zero; A(0,2) = zero;
235 A(1,0) = zero; A(1,1) = gamma; A(1,2) = zero;
236 A(2,0) = zero; A(2,1) = (1.0-2.0*gamma); A(2,2) = gamma;
239 b(0) = zero; b(1) = onehalf; b(2) = onehalf;
242 c(0) = zero; c(1) = gamma; c(2) = one-gamma;
247 "Implicit Tableau - Partitioned IMEX RK ARS 233",
248 A,b,c,order,order,order));
250 this->setImplicitTableau(impTableau);
252 this->setStepperType(
"Partitioned IMEX RK ARS 233");
255 }
else if (stepperType ==
"General Partitioned IMEX RK") {
256 this->setExplicitTableau(explicitTableau);
257 this->setImplicitTableau(implicitTableau);
258 this->setStepperType(
"General Partitioned IMEX RK");
262 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
263 "Error - Not a valid StepperIMEX_RK_Partition type! Stepper Type = "
264 << stepperType <<
"\n"
265 <<
" Current valid types are: " <<
"\n"
266 <<
" 'Partitioned IMEX RK 1st order'" <<
"\n"
267 <<
" 'Partitioned IMEX RK SSP2'" <<
"\n"
268 <<
" 'Partitioned IMEX RK ARS 233'" <<
"\n"
269 <<
" 'General Partitioned IMEX RK'" <<
"\n");
272 TEUCHOS_TEST_FOR_EXCEPTION(explicitTableau_==Teuchos::null,
274 "Error - StepperIMEX_RK_Partition - Explicit tableau is null!");
275 TEUCHOS_TEST_FOR_EXCEPTION(implicitTableau_==Teuchos::null,
277 "Error - StepperIMEX_RK_Partition - Implicit tableau is null!");
278 TEUCHOS_TEST_FOR_EXCEPTION(
279 explicitTableau_->numStages()!=implicitTableau_->numStages(),
281 "Error - StepperIMEX_RK_Partition - Number of stages do not match!\n"
282 <<
" Explicit tableau = " << explicitTableau_->description() <<
"\n"
283 <<
" number of stages = " << explicitTableau_->numStages() <<
"\n"
284 <<
" Implicit tableau = " << implicitTableau_->description() <<
"\n"
285 <<
" number of stages = " << implicitTableau_->numStages() <<
"\n");
287 this->isInitialized_ =
false;
291 template<
class Scalar>
295 TEUCHOS_TEST_FOR_EXCEPTION(explicitTableau->isImplicit() ==
true,
297 "Error - Received an implicit Tableau for setExplicitTableau()!\n" <<
298 " Tableau = " << explicitTableau->description() <<
"\n");
299 explicitTableau_ = explicitTableau;
301 this->isInitialized_ =
false;
305 template<
class Scalar>
309 TEUCHOS_TEST_FOR_EXCEPTION( implicitTableau->isDIRK() !=
true,
311 "Error - Did not receive a DIRK Tableau for setImplicitTableau()!\n" <<
312 " Tableau = " << implicitTableau->description() <<
"\n");
313 implicitTableau_ = implicitTableau;
315 this->isInitialized_ =
false;
318 template<
class Scalar>
320 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& appModel)
323 using Teuchos::rcp_const_cast;
324 using Teuchos::rcp_dynamic_cast;
325 RCP<Thyra::ModelEvaluator<Scalar> > ncModel =
326 rcp_const_cast<Thyra::ModelEvaluator<Scalar> > (appModel);
327 RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> > modelPairIMEX =
329 TEUCHOS_TEST_FOR_EXCEPTION( modelPairIMEX == Teuchos::null, std::logic_error,
330 "Error - StepperIMEX_RK::setModel() was given a ModelEvaluator that\n"
331 " could not be cast to a WrapperModelEvaluatorPairPartIMEX_Basic!\n"
332 " From: " << appModel <<
"\n"
333 " To : " << modelPairIMEX <<
"\n"
334 " Likely have given the wrong ModelEvaluator to this Stepper.\n");
336 setModelPair(modelPairIMEX);
338 this->isInitialized_ =
false;
346 template<
class Scalar>
351 Teuchos::RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> >
352 wrapperModelPairIMEX =
354 (this->wrapperModel_);
357 wrapperModelPairIMEX = mePairIMEX;
358 wrapperModelPairIMEX->initialize();
360 this->wrapperModel_ = wrapperModelPairIMEX;
362 this->isInitialized_ =
false;
370 template<
class Scalar>
372 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& explicitModel,
373 const Teuchos::RCP<
const Thyra::ModelEvaluator<Scalar> >& implicitModel)
377 this->wrapperModel_ = Teuchos::rcp(
379 explicitModel, implicitModel));
381 this->isInitialized_ =
false;
385 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
386 template<
class Scalar>
390 if (this->stepperObserver_ == Teuchos::null)
391 this->stepperObserver_ =
394 if (( obs == Teuchos::null ) and (this->stepperObserver_->getSize() >0 ) )
397 if (( obs == Teuchos::null ) and (this->stepperObserver_->getSize() == 0) )
402 this->stepperObserver_->addObserver(
405 Teuchos::RCP<Teuchos::FancyOStream> out = this->getOStream();
406 Teuchos::OSTab ostab(out,0,
"setObserver");
407 *out <<
"Tempus::StepperIMEX_RK_Partition::setObserver: Warning: An observer has been provided that";
408 *out <<
" does not support Tempus::StepperRKObserver. This observer WILL NOT be added.";
409 *out <<
" In the future, this will result in a runtime error!" << std::endl;
412 this->isInitialized_ =
false;
417 template<
class Scalar>
420 Teuchos::RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> >
421 wrapperModelPairIMEX =
423 (this->wrapperModel_);
424 TEUCHOS_TEST_FOR_EXCEPTION( wrapperModelPairIMEX == Teuchos::null,
426 "Error - Can not cast the wrapper Model Evaluator to a IMEX Model Pair."
427 "StepperIMEX_RK_Partition::initialize()\n");
430 const int numStages = explicitTableau_->numStages();
431 stageF_.resize(numStages);
432 stageGx_.resize(numStages);
433 for(
int i=0; i < numStages; i++) {
434 stageF_[i] = Thyra::createMember(wrapperModelPairIMEX->
435 getExplicitModel()->get_f_space());
436 stageGx_[i] = Thyra::createMember(wrapperModelPairIMEX->
437 getImplicitModel()->get_f_space());
438 assign(stageF_[i].ptr(), Teuchos::ScalarTraits<Scalar>::zero());
439 assign(stageGx_[i].ptr(), Teuchos::ScalarTraits<Scalar>::zero());
442 xTilde_ = Thyra::createMember(wrapperModelPairIMEX->
443 getImplicitModel()->get_x_space());
444 assign(xTilde_.ptr(), Teuchos::ScalarTraits<Scalar>::zero());
450 template<
class Scalar>
456 int numStates = solutionHistory->getNumStates();
458 TEUCHOS_TEST_FOR_EXCEPTION(numStates < 1, std::logic_error,
459 "Error - setInitialConditions() needs at least one SolutionState\n"
460 " to set the initial condition. Number of States = " << numStates);
463 RCP<Teuchos::FancyOStream> out = this->getOStream();
464 Teuchos::OSTab ostab(out,1,
"StepperIMEX_RK::setInitialConditions()");
465 *out <<
"Warning -- SolutionHistory has more than one state!\n"
466 <<
"Setting the initial conditions on the currentState.\n"<<std::endl;
469 RCP<SolutionState<Scalar> > initialState = solutionHistory->getCurrentState();
470 RCP<Thyra::VectorBase<Scalar> > x = initialState->getX();
473 auto inArgs = this->wrapperModel_->getNominalValues();
474 if (x == Teuchos::null) {
475 TEUCHOS_TEST_FOR_EXCEPTION( (x == Teuchos::null) &&
476 (inArgs.get_x() == Teuchos::null), std::logic_error,
477 "Error - setInitialConditions() needs the ICs from the SolutionHistory\n"
478 " or getNominalValues()!\n");
480 x = Teuchos::rcp_const_cast<Thyra::VectorBase<Scalar> >(inArgs.get_x());
481 initialState->setX(x);
486 std::string icConsistency = this->getICConsistency();
487 TEUCHOS_TEST_FOR_EXCEPTION(icConsistency !=
"None", std::logic_error,
488 "Error - setInitialConditions() requested a consistency of '"
489 << icConsistency <<
"'.\n"
490 " But only 'None' is available for IMEX-RK!\n");
492 TEUCHOS_TEST_FOR_EXCEPTION( this->getUseFSAL(), std::logic_error,
493 "Error - The First-Step-As-Last (FSAL) principle is not "
494 <<
"available for IMEX-RK. Set useFSAL=false.\n");
498 template <
typename Scalar>
500 const Teuchos::RCP<
const Thyra::VectorBase<Scalar> > & X,
501 const Teuchos::RCP<
const Thyra::VectorBase<Scalar> > & Y,
502 Scalar time, Scalar stepSize, Scalar stageNumber,
503 const Teuchos::RCP<Thyra::VectorBase<Scalar> > & G)
const
505 typedef Thyra::ModelEvaluatorBase MEB;
506 Teuchos::RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> >
507 wrapperModelPairIMEX =
509 (this->wrapperModel_);
510 MEB::InArgs<Scalar> inArgs = wrapperModelPairIMEX->
getInArgs();
512 inArgs.set_p(wrapperModelPairIMEX->getParameterIndex(), Y);
513 if (inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(time);
514 if (inArgs.supports(MEB::IN_ARG_step_size)) inArgs.set_step_size(stepSize);
515 if (inArgs.supports(MEB::IN_ARG_stage_number))
516 inArgs.set_stage_number(stageNumber);
523 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(Teuchos::null);
525 MEB::OutArgs<Scalar> outArgs = wrapperModelPairIMEX->getOutArgs();
528 wrapperModelPairIMEX->getImplicitModel()->evalModel(inArgs,outArgs);
529 Thyra::Vt_S(G.ptr(), -1.0);
533 template <
typename Scalar>
535 const Teuchos::RCP<
const Thyra::VectorBase<Scalar> > & Z,
536 Scalar time, Scalar stepSize, Scalar stageNumber,
537 const Teuchos::RCP<Thyra::VectorBase<Scalar> > & F)
const
539 typedef Thyra::ModelEvaluatorBase MEB;
541 Teuchos::RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> >
542 wrapperModelPairIMEX =
544 (this->wrapperModel_);
545 MEB::InArgs<Scalar> inArgs =
548 if (inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(time);
549 if (inArgs.supports(MEB::IN_ARG_step_size)) inArgs.set_step_size(stepSize);
550 if (inArgs.supports(MEB::IN_ARG_stage_number))
551 inArgs.set_stage_number(stageNumber);
558 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(Teuchos::null);
560 MEB::OutArgs<Scalar> outArgs =
561 wrapperModelPairIMEX->getExplicitModel()->createOutArgs();
564 wrapperModelPairIMEX->getExplicitModel()->evalModel(inArgs, outArgs);
565 Thyra::Vt_S(F.ptr(), -1.0);
569 template<
class Scalar>
573 this->checkInitialized();
576 using Teuchos::SerialDenseMatrix;
577 using Teuchos::SerialDenseVector;
579 TEMPUS_FUNC_TIME_MONITOR(
"Tempus::StepperIMEX_RK_Partition::takeStep()");
581 TEUCHOS_TEST_FOR_EXCEPTION(solutionHistory->getNumStates() < 2,
583 "Error - StepperIMEX_RK_Partition<Scalar>::takeStep(...)\n"
584 "Need at least two SolutionStates for IMEX_RK_Partition.\n"
585 " Number of States = " << solutionHistory->getNumStates() <<
"\n"
586 "Try setting in \"Solution History\" \"Storage Type\" = \"Undo\"\n"
587 " or \"Storage Type\" = \"Static\" and \"Storage Limit\" = \"2\"\n");
589 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
590 this->stepperObserver_->observeBeginTakeStep(solutionHistory, *
this);
592 RCP<StepperIMEX_RK_Partition<Scalar> > thisStepper = Teuchos::rcpFromRef(*
this);
593 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
596 RCP<SolutionState<Scalar> > currentState=solutionHistory->getCurrentState();
597 RCP<SolutionState<Scalar> > workingState=solutionHistory->getWorkingState();
598 const Scalar dt = workingState->getTimeStep();
599 const Scalar time = currentState->getTime();
601 const int numStages = explicitTableau_->numStages();
602 const SerialDenseMatrix<int,Scalar> & AHat = explicitTableau_->A();
603 const SerialDenseVector<int,Scalar> & bHat = explicitTableau_->b();
604 const SerialDenseVector<int,Scalar> & cHat = explicitTableau_->c();
605 const SerialDenseMatrix<int,Scalar> & A = implicitTableau_->A();
606 const SerialDenseVector<int,Scalar> & b = implicitTableau_->b();
607 const SerialDenseVector<int,Scalar> & c = implicitTableau_->c();
609 Teuchos::RCP<WrapperModelEvaluatorPairPartIMEX_Basic<Scalar> >
610 wrapperModelPairIMEX =
612 (this->wrapperModel_);
615 Thyra::SolveStatus<Scalar> sStatus;
616 stageZ_ = workingState->getX();
617 Thyra::assign(stageZ_.ptr(), *(currentState->getX()));
618 RCP<Thyra::VectorBase<Scalar> > stageY =
619 wrapperModelPairIMEX->getExplicitOnlyVector(stageZ_);
620 RCP<Thyra::VectorBase<Scalar> > stageX =
621 wrapperModelPairIMEX->getIMEXVector(stageZ_);
624 for (
int i = 0; i < numStages; ++i) {
625 this->setStageNumber(i);
626 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
627 this->stepperObserver_->observeBeginStage(solutionHistory, *
this);
630 Thyra::assign(stageY.ptr(),
631 *(wrapperModelPairIMEX->getExplicitOnlyVector(currentState->getX())));
632 Thyra::assign(xTilde_.ptr(),
633 *(wrapperModelPairIMEX->getIMEXVector(currentState->getX())));
634 for (
int j = 0; j < i; ++j) {
635 if (AHat(i,j) != Teuchos::ScalarTraits<Scalar>::zero()) {
636 RCP<Thyra::VectorBase<Scalar> > stageFy =
637 wrapperModelPairIMEX->getExplicitOnlyVector(stageF_[j]);
638 RCP<Thyra::VectorBase<Scalar> > stageFx =
639 wrapperModelPairIMEX->getIMEXVector(stageF_[j]);
640 Thyra::Vp_StV(stageY.ptr(), -dt*AHat(i,j), *stageFy);
641 Thyra::Vp_StV(xTilde_.ptr(), -dt*AHat(i,j), *stageFx);
643 if (A (i,j) != Teuchos::ScalarTraits<Scalar>::zero())
644 Thyra::Vp_StV(xTilde_.ptr(), -dt*A (i,j), *(stageGx_[j]));
647 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
650 Scalar ts = time + c(i)*dt;
651 Scalar tHats = time + cHat(i)*dt;
652 if (A(i,i) == Teuchos::ScalarTraits<Scalar>::zero()) {
654 bool isNeeded =
false;
655 for (
int k=i+1; k<numStages; ++k)
if (A(k,i) != 0.0) isNeeded =
true;
656 if (b(i) != 0.0) isNeeded =
true;
657 if (isNeeded ==
false) {
659 assign(stageGx_[i].ptr(), Teuchos::ScalarTraits<Scalar>::zero());
661 Thyra::assign(stageX.ptr(), *xTilde_);
662 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
663 this->stepperObserver_->observeBeforeImplicitExplicitly(solutionHistory, *
this);
665 evalImplicitModelExplicitly(stageX, stageY, ts, dt, i, stageGx_[i]);
669 const Scalar alpha = Scalar(1.0)/(dt*A(i,i));
670 const Scalar beta = Scalar(1.0);
673 RCP<TimeDerivative<Scalar> > timeDer =
675 alpha, xTilde_.getConst()));
678 typedef Thyra::ModelEvaluatorBase MEB;
681 wrapperModelPairIMEX->setUseImplicitModel(
true);
682 MEB::InArgs<Scalar> inArgs = wrapperModelPairIMEX->createInArgs();
683 MEB::OutArgs<Scalar> outArgs = wrapperModelPairIMEX->createOutArgs();
684 inArgs.set_x(stageX);
685 if (wrapperModelPairIMEX->getParameterIndex() >= 0)
686 inArgs.set_p(wrapperModelPairIMEX->getParameterIndex(), stageY);
687 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(stageGx_[i]);
688 if (inArgs.supports(MEB::IN_ARG_t )) inArgs.set_t (ts);
689 if (inArgs.supports(MEB::IN_ARG_step_size)) inArgs.set_step_size(dt);
690 if (inArgs.supports(MEB::IN_ARG_alpha )) inArgs.set_alpha (alpha);
691 if (inArgs.supports(MEB::IN_ARG_beta )) inArgs.set_beta (beta);
692 if (inArgs.supports(MEB::IN_ARG_stage_number))
693 inArgs.set_stage_number(i);
695 wrapperModelPairIMEX->setForSolve(timeDer, inArgs, outArgs);
697 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
698 this->stepperObserver_->observeBeforeSolve(solutionHistory, *
this);
700 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
703 this->solver_->setModel(wrapperModelPairIMEX);
704 sStatus = this->solveImplicitODE(stageX);
705 if (sStatus.solveStatus != Thyra::SOLVE_STATUS_CONVERGED) pass =
false;
707 wrapperModelPairIMEX->setUseImplicitModel(
false);
709 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
710 this->stepperObserver_->observeAfterSolve(solutionHistory, *
this);
712 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
716 Thyra::V_StVpStV(stageGx_[i].ptr(), -alpha, *stageX, alpha, *xTilde_);
719 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
720 this->stepperObserver_->observeBeforeExplicit(solutionHistory, *
this);
722 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
724 evalExplicitModel(stageZ_, tHats, dt, i, stageF_[i]);
725 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
726 this->stepperObserver_->observeEndStage(solutionHistory, *
this);
728 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
734 Thyra::assign((workingState->getX()).ptr(), *(currentState->getX()));
735 RCP<Thyra::VectorBase<Scalar> > Z = workingState->getX();
736 RCP<Thyra::VectorBase<Scalar> > X = wrapperModelPairIMEX->getIMEXVector(Z);
737 for (
int i=0; i < numStages; ++i) {
738 if (bHat(i) != Teuchos::ScalarTraits<Scalar>::zero())
739 Thyra::Vp_StV(Z.ptr(), -dt*bHat(i), *(stageF_[i]));
740 if (b (i) != Teuchos::ScalarTraits<Scalar>::zero())
741 Thyra::Vp_StV(X.ptr(), -dt*b (i), *(stageGx_[i]));
744 if (pass ==
true) workingState->setSolutionStatus(
Status::PASSED);
746 workingState->setOrder(this->getOrder());
747 workingState->computeNorms(currentState);
748 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
749 this->stepperObserver_->observeEndTakeStep(solutionHistory, *
this);
751 this->stepperRKAppAction_->execute(solutionHistory, thisStepper,
755 this->setStageNumber(-1);
765 template<
class Scalar>
766 Teuchos::RCP<Tempus::StepperState<Scalar> >
770 Teuchos::RCP<Tempus::StepperState<Scalar> > stepperState =
776 template<
class Scalar>
778 Teuchos::FancyOStream &out,
779 const Teuchos::EVerbosityLevel verbLevel)
const
785 out <<
"--- StepperIMEX_RK_Partition ---\n";
786 out <<
" explicitTableau_ = " << explicitTableau_ << std::endl;
787 if (verbLevel == Teuchos::VERB_HIGH)
788 explicitTableau_->describe(out, verbLevel);
789 out <<
" implicitTableau_ = " << implicitTableau_ << std::endl;
790 if (verbLevel == Teuchos::VERB_HIGH)
791 implicitTableau_->describe(out, verbLevel);
792 out <<
" xTilde_ = " << xTilde_ << std::endl;
793 out <<
" stageZ_ = " << stageZ_ << std::endl;
794 out <<
" stageF_.size() = " << stageF_.size() << std::endl;
795 int numStages = stageF_.size();
796 for (
int i=0; i<numStages; ++i)
797 out <<
" stageF_["<<i<<
"] = " << stageF_[i] << std::endl;
798 out <<
" stageGx_.size() = " << stageGx_.size() << std::endl;
799 numStages = stageGx_.size();
800 for (
int i=0; i<numStages; ++i)
801 out <<
" stageGx_["<<i<<
"] = " << stageGx_[i] << std::endl;
802 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
803 out <<
" stepperObserver_ = " << stepperObserver_ << std::endl;
805 out <<
" stepperRKAppAction_= " << this->stepperRKAppAction_ << std::endl;
806 out <<
" order_ = " << order_ << std::endl;
807 out <<
"--------------------------------" << std::endl;
811 template<
class Scalar>
814 bool isValidSetup =
true;
818 Teuchos::RCP<WrapperModelEvaluatorPairIMEX<Scalar> > wrapperModelPairIMEX =
820 this->wrapperModel_);
822 if ( wrapperModelPairIMEX->getExplicitModel() == Teuchos::null) {
823 isValidSetup =
false;
824 out <<
"The explicit ModelEvaluator is not set!\n";
827 if ( wrapperModelPairIMEX->getImplicitModel() == Teuchos::null) {
828 isValidSetup =
false;
829 out <<
"The implicit ModelEvaluator is not set!\n";
832 if (this->wrapperModel_ == Teuchos::null) {
833 isValidSetup =
false;
834 out <<
"The wrapper ModelEvaluator is not set!\n";
837 if (this->solver_ == Teuchos::null) {
838 isValidSetup =
false;
839 out <<
"The solver is not set!\n";
842 #ifndef TEMPUS_HIDE_DEPRECATED_CODE
843 if (stepperObserver_ == Teuchos::null) {
844 isValidSetup =
false;
845 out <<
"The stepper observer is not set!\n";
848 if (this->stepperRKAppAction_ == Teuchos::null) {
849 isValidSetup =
false;
850 out <<
"The AppAction is not set!\n";
853 if ( explicitTableau_ == Teuchos::null ) {
854 isValidSetup =
false;
855 out <<
"The explicit tableau is not set!\n";
858 if ( implicitTableau_ == Teuchos::null ) {
859 isValidSetup =
false;
860 out <<
"The implicit tableau is not set!\n";
867 template<
class Scalar>
868 Teuchos::RCP<const Teuchos::ParameterList>
871 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
872 pl->setName(
"Default Stepper - Partitioned IMEX RK SSP2");
873 pl->set<std::string>(
"Stepper Type",
"Partitioned IMEX RK SSP2");
875 pl->set<
bool>(
"Initial Condition Consistency Check",
876 this->getICConsistencyCheckDefault());
877 pl->set<std::string>(
"Solver Name",
"Default Solver");
878 pl->set<
bool> (
"Zero Initial Guess",
false);
880 pl->set(
"Default Solver", *solverPL);
887 #endif // Tempus_StepperIMEX_RK_Partition_impl_hpp
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
StepperIMEX_RK_Partition()
Default constructor.
virtual void setTableaus(std::string stepperType="", Teuchos::RCP< const RKButcherTableau< Scalar > > explicitTableau=Teuchos::null, Teuchos::RCP< const RKButcherTableau< Scalar > > implicitTableau=Teuchos::null)
Set both the explicit and implicit tableau from ParameterList.
virtual void setModelPair(const Teuchos::RCP< WrapperModelEvaluatorPairPartIMEX_Basic< Scalar > > &modelPair)
Create WrapperModelPairIMEX from user-supplied ModelEvaluator pair.
ModelEvaluator pair for implicit and explicit (IMEX) evaluations.
Teuchos::RCP< Teuchos::ParameterList > defaultSolverParameters()
Returns the default solver ParameterList for implicit Steppers.
void validExplicitODE(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Validate that the model supports explicit ODE evaluation, f(x,t) [=xdot].
virtual void initialize()
Initialize after construction and changing input parameters.
virtual Thyra::ModelEvaluatorBase::InArgs< Scalar > getInArgs()
Get InArgs the wrapper ModelEvalutor.
void evalExplicitModel(const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &X, Scalar time, Scalar stepSize, Scalar stageNumber, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &F) const
Thyra Base interface for time steppers.
StepperState is a simple class to hold state information about the stepper.
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Provide a StepperState to the SolutionState. This Stepper does not have any special state data...
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Application Action for StepperRKBase.
void evalImplicitModelExplicitly(const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &X, const Teuchos::RCP< const Thyra::VectorBase< Scalar > > &Y, Scalar time, Scalar stepSize, Scalar stageNumber, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &G) const
StepperObserver class for Stepper class.
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual void setExplicitTableau(Teuchos::RCP< const RKButcherTableau< Scalar > > explicitTableau)
Set the explicit tableau from tableau.
virtual void setObserver(Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null)
Set Observer.
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
virtual void setImplicitTableau(Teuchos::RCP< const RKButcherTableau< Scalar > > implicitTableau)
Set the implicit tableau from tableau.
Time-derivative interface for Partitioned IMEX RK.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
This observer is a composite observer,.
virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > getExplicitModel() const
StepperRKObserver class for StepperRK.
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
void validImplicitODE_DAE(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Validate ME supports implicit ODE/DAE evaluation, f(xdot,x,t) [= 0].
virtual void initialize()
Initialize during construction and after changing input parameters.
void getValidParametersBasic(Teuchos::RCP< Teuchos::ParameterList > pl, std::string stepperType)
Provide basic parameters to Steppers.