29 #ifndef RYTHMOS_DEFAULT_INTEGRATOR_DEF_HPP 
   30 #define RYTHMOS_DEFAULT_INTEGRATOR_DEF_HPP 
   32 #include "Rythmos_DefaultIntegrator_decl.hpp" 
   33 #include "Rythmos_InterpolationBufferHelpers.hpp" 
   34 #include "Rythmos_IntegrationControlStrategyBase.hpp" 
   35 #include "Rythmos_IntegrationObserverBase.hpp" 
   36 #include "Rythmos_InterpolationBufferAppenderBase.hpp" 
   37 #include "Rythmos_PointwiseInterpolationBufferAppender.hpp" 
   38 #include "Rythmos_StepperHelpers.hpp" 
   39 #include "Teuchos_VerboseObjectParameterListHelpers.hpp" 
   40 #include "Teuchos_Assert.hpp" 
   41 #include "Teuchos_as.hpp" 
   50 template<
class Scalar>
 
   51 RCP<DefaultIntegrator<Scalar> >
 
   54   RCP<DefaultIntegrator<Scalar> >
 
   64 template<
class Scalar>
 
   65 RCP<DefaultIntegrator<Scalar> >
 
   71   RCP<DefaultIntegrator<Scalar> >
 
   73   integrator->setIntegrationControlStrategy(integrationControlStrategy);
 
   74   integrator->setIntegrationObserver(integrationObserver);
 
   83 template<
class Scalar>
 
   84 RCP<DefaultIntegrator<Scalar> >
 
   89   RCP<DefaultIntegrator<Scalar> >
 
   91   integrator->setIntegrationControlStrategy(integrationControlStrategy);
 
  100 template<
class Scalar>
 
  101 RCP<DefaultIntegrator<Scalar> >
 
  106   RCP<DefaultIntegrator<Scalar> >
 
  108   integrator->setIntegrationObserver(integrationObserver);
 
  121 template<
class Scalar>
 
  125 template<
class Scalar>
 
  128   std::numeric_limits<int>::max();
 
  135 template<
class Scalar>
 
  137   :landOnFinalTime_(true),
 
  138    maxNumTimeSteps_(maxNumTimeSteps_default_),
 
  139    currTimeStepIndex_(-1)
 
  143 template<
class Scalar>
 
  148 #ifdef HAVE_RYTHMOS_DEBUG 
  149   TEUCHOS_TEST_FOR_EXCEPT(is_null(integrationControlStrategy));
 
  151   integrationControlStrategy_ = integrationControlStrategy;
 
  154 template<
class Scalar>
 
  155 RCP<IntegrationControlStrategyBase<Scalar> >
 
  158   return integrationControlStrategy_;
 
  161 template<
class Scalar>
 
  162 RCP<const IntegrationControlStrategyBase<Scalar> >
 
  165   return integrationControlStrategy_;
 
  169 template<
class Scalar>
 
  174 #ifdef HAVE_RYTHMOS_DEBUG 
  175   TEUCHOS_TEST_FOR_EXCEPT(is_null(integrationObserver));
 
  177   integrationObserver_ = integrationObserver;
 
  181 template<
class Scalar>
 
  186   interpBufferAppender_ = interpBufferAppender.assert_not_null();
 
  190 template<
class Scalar>
 
  191 RCP<const InterpolationBufferAppenderBase<Scalar> >
 
  194   return interpBufferAppender_;
 
  197 template<
class Scalar>
 
  198 RCP<InterpolationBufferAppenderBase<Scalar> >
 
  201   return interpBufferAppender_;
 
  204 template<
class Scalar>
 
  205 RCP<InterpolationBufferAppenderBase<Scalar> >
 
  208   RCP<InterpolationBufferAppenderBase<Scalar> > InterpBufferAppender;
 
  209   std::swap( InterpBufferAppender, interpBufferAppender_ );
 
  210   return InterpBufferAppender;
 
  217 template<
class Scalar>
 
  219   RCP<ParameterList> 
const& paramList
 
  222   TEUCHOS_TEST_FOR_EXCEPT(is_null(paramList));
 
  223   paramList->validateParameters(*getValidParameters());
 
  224   this->setMyParamList(paramList);
 
  225   maxNumTimeSteps_ = paramList->get(
 
  226     maxNumTimeSteps_name_, maxNumTimeSteps_default_);
 
  227   Teuchos::readVerboseObjectSublist(&*paramList,
this);
 
  231 template<
class Scalar>
 
  232 RCP<const ParameterList>
 
  235   static RCP<const ParameterList> validPL;
 
  236   if (is_null(validPL) ) {
 
  237     RCP<ParameterList> pl = Teuchos::parameterList();
 
  238     pl->set(maxNumTimeSteps_name_, maxNumTimeSteps_default_,
 
  239       "Set the maximum number of integration time-steps allowed.");
 
  240     Teuchos::setupVerboseObjectSublist(&*pl);
 
  250 template<
class Scalar>
 
  251 RCP<IntegratorBase<Scalar> >
 
  256   RCP<DefaultIntegrator<Scalar> >
 
  259   newIntegrator->stepper_ = Teuchos::null;
 
  260   const RCP<const ParameterList> paramList = this->getParameterList();
 
  261   if (!is_null(paramList)) {
 
  262     newIntegrator->setParameterList(Teuchos::parameterList(*paramList));
 
  264   if (!is_null(integrationControlStrategy_)) {
 
  265     newIntegrator->integrationControlStrategy_ =
 
  266       integrationControlStrategy_->cloneIntegrationControlStrategy().assert_not_null();
 
  268   if (!is_null(integrationObserver_)) {
 
  269     newIntegrator->integrationObserver_ =
 
  270       integrationObserver_->cloneIntegrationObserver().assert_not_null();
 
  272   if (!is_null(trailingInterpBuffer_)) {
 
  274     newIntegrator->trailingInterpBuffer_ = null;
 
  278   if (!is_null(interpBufferAppender_)) {
 
  280     newIntegrator->interpBufferAppender_ = null;
 
  284   return newIntegrator;
 
  288 template<
class Scalar>
 
  291   const Scalar &finalTime,
 
  292   const bool landOnFinalTime
 
  295   typedef Teuchos::ScalarTraits<Scalar> ST;
 
  296   TEUCHOS_TEST_FOR_EXCEPT(is_null(stepper));
 
  297   TEUCHOS_TEST_FOR_EXCEPT( finalTime < stepper->getTimeRange().lower() );
 
  298   TEUCHOS_ASSERT( stepper->getTimeRange().length() == ST::zero() );
 
  301   integrationTimeDomain_ = timeRange(stepper_->getTimeRange().lower(),
 
  303   landOnFinalTime_ = landOnFinalTime;
 
  304   currTimeStepIndex_ = 0;
 
  306   if (!is_null(integrationControlStrategy_))
 
  307     integrationControlStrategy_->resetIntegrationControlStrategy(
 
  308       integrationTimeDomain_
 
  310   if (!is_null(integrationObserver_))
 
  311     integrationObserver_->resetIntegrationObserver(
 
  312       integrationTimeDomain_
 
  317 template<
class Scalar>
 
  320   RCP<StepperBase<Scalar> > stepper_temp = stepper_;
 
  321   stepper_ = Teuchos::null;
 
  322   return(stepper_temp);
 
  326 template<
class Scalar>
 
  333 template<
class Scalar>
 
  340 template<
class Scalar>
 
  345   trailingInterpBuffer_ = trailingInterpBuffer.assert_not_null();
 
  349 template<
class Scalar>
 
  350 RCP<InterpolationBufferBase<Scalar> >
 
  353   return trailingInterpBuffer_;
 
  357 template<
class Scalar>
 
  358 RCP<const InterpolationBufferBase<Scalar> >
 
  361   return trailingInterpBuffer_;
 
  364 template<
class Scalar>
 
  365 RCP<InterpolationBufferBase<Scalar> >
 
  368   RCP<InterpolationBufferBase<Scalar> > trailingInterpBuffer;
 
  369   std::swap( trailingInterpBuffer, trailingInterpBuffer_ );
 
  370   return trailingInterpBuffer;
 
  374 template<
class Scalar>
 
  376   const Array<Scalar>& time_vec,
 
  377   Array<RCP<
const Thyra::VectorBase<Scalar> > >* x_vec,
 
  378   Array<RCP<
const Thyra::VectorBase<Scalar> > >* xdot_vec,
 
  379   Array<ScalarMag>* accuracy_vec)
 
  382   RYTHMOS_FUNC_TIME_MONITOR_DIFF(
"Rythmos:DefaultIntegrator::getFwdPoints",
 
  385   using Teuchos::incrVerbLevel;
 
  387   using Teuchos::Describable;
 
  391   typedef Teuchos::VerboseObjectTempState<IBB> VOTSIBB;
 
  395   RCP<Teuchos::FancyOStream> out = this->getOStream();
 
  396   Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
 
  397   Teuchos::OSTab tab(out);
 
  398   VOTSIBB stepper_outputTempState(stepper_,out,incrVerbLevel(verbLevel,-1));
 
  400   if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  401     *out << 
"\nEntering " << this->Describable::description()
 
  402          << 
"::getFwdPoints(...) ...\n" 
  403          << 
"\nStepper: " << Teuchos::describe(*stepper_,verbLevel);
 
  405   if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
 
  406     *out << 
"\nRequested time points: " << Teuchos::toString(time_vec) << 
"\n";
 
  409   if (!is_null(integrationObserver_)) {
 
  410     integrationObserver_->setOStream(out);
 
  411     integrationObserver_->setVerbLevel(incrVerbLevel(verbLevel,-1));
 
  412     integrationObserver_->observeStartTimeIntegration(*stepper_);
 
  419   const int numTimePoints = time_vec.size();
 
  422   assertTimePointsAreSorted(time_vec);
 
  423   TEUCHOS_TEST_FOR_EXCEPT(accuracy_vec!=0); 
 
  427     x_vec->resize(numTimePoints);
 
  429     xdot_vec->resize(numTimePoints);
 
  434   int nextTimePointIndex = 0;
 
  436   assertNoTimePointsBeforeCurrentTimeRange(*
this,time_vec,nextTimePointIndex);
 
  443     RYTHMOS_FUNC_TIME_MONITOR(
 
  444       "Rythmos:DefaultIntegrator::getFwdPoints: getPoints");
 
  446     getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex);
 
  454   while ( nextTimePointIndex < numTimePoints ) {
 
  458     const Scalar t = time_vec[nextTimePointIndex];
 
  459     bool advanceStepperToTimeSucceeded = 
false;
 
  462     if ( integrationTimeDomain_.upper() < t ) {
 
  463       *out << 
"\nWARNING: Skipping time point, t = " << t
 
  464            << 
", because it is beyond 'Final Time' = " 
  465            << integrationTimeDomain_.upper() << 
"\n";
 
  469       RYTHMOS_FUNC_TIME_MONITOR(
 
  470         "Rythmos:DefaultIntegrator::getFwdPoints: advanceStepperToTime");
 
  471       advanceStepperToTimeSucceeded = advanceStepperToTime(t);
 
  473     if (!advanceStepperToTimeSucceeded) {
 
  474       bool reachedMaxNumTimeSteps = (currTimeStepIndex_ >= maxNumTimeSteps_);
 
  475       if (reachedMaxNumTimeSteps) {
 
  479       TEUCHOS_TEST_FOR_EXCEPTION(
 
  481           this->description() << 
"\n\n" 
  482           "Error:  The integration failed to get to time " << t
 
  483           << 
" and only achieved\ngetting to " 
  484           << stepper_->getTimeRange().upper() << 
"!" 
  490       RYTHMOS_FUNC_TIME_MONITOR(
 
  491         "Rythmos:DefaultIntegrator::getFwdPoints: getPoints (fwd)");
 
  492       getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex);
 
  498   if (!is_null(integrationObserver_)) {
 
  499     integrationObserver_->observeEndTimeIntegration(*stepper_);
 
  502   if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  503     *out << 
"\nLeaving " << this->Describable::description()
 
  504          << 
"::getFwdPoints(...) ...\n";
 
  508 template<
class Scalar>
 
  513     stepper_->getTimeRange().lower(),
 
  514     integrationTimeDomain_.upper()
 
  522 template<
class Scalar>
 
  523 RCP<const Thyra::VectorSpaceBase<Scalar> >
 
  526   return stepper_->get_x_space();
 
  530 template<
class Scalar>
 
  532   const Array<Scalar>& time_vec,
 
  533   const Array<RCP<
const Thyra::VectorBase<Scalar> > >& x_vec,
 
  534   const Array<RCP<
const Thyra::VectorBase<Scalar> > >& xdot_vec
 
  537   stepper_->addPoints(time_vec,x_vec,xdot_vec);
 
  541 template<
class Scalar>
 
  543   const Array<Scalar>& time_vec,
 
  544   Array<RCP<
const Thyra::VectorBase<Scalar> > >* x_vec,
 
  545   Array<RCP<
const Thyra::VectorBase<Scalar> > >* xdot_vec,
 
  546   Array<ScalarMag>* accuracy_vec
 
  561   stepper_->getPoints(time_vec, x_vec, xdot_vec, accuracy_vec);
 
  566 template<
class Scalar>
 
  569   if (nonnull(trailingInterpBuffer_)) {
 
  570     return timeRange(trailingInterpBuffer_->getTimeRange().lower(),
 
  571       stepper_->getTimeRange().upper());
 
  573   return stepper_->getTimeRange();
 
  577 template<
class Scalar>
 
  580   stepper_->getNodes(time_vec);
 
  584 template<
class Scalar>
 
  587   stepper_->removeNodes(time_vec);
 
  591 template<
class Scalar>
 
  594   return stepper_->getOrder();
 
  601 template<
class Scalar>
 
  604   if (!is_null(trailingInterpBuffer_) && is_null(interpBufferAppender_))
 
  605     interpBufferAppender_ = pointwiseInterpolationBufferAppender<Scalar>();
 
  610 template<
class Scalar>
 
  611 bool DefaultIntegrator<Scalar>::advanceStepperToTime(
const Scalar& advance_to_t)
 
  614   RYTHMOS_FUNC_TIME_MONITOR_DIFF(
 
  615     "Rythmos:DefaultIntegrator::advanceStepperToTime", TopLevel);
 
  618   typedef std::numeric_limits<Scalar> NL;
 
  619   using Teuchos::incrVerbLevel;
 
  621   using Teuchos::Describable;
 
  623   using Teuchos::OSTab;
 
  624   typedef Teuchos::ScalarTraits<Scalar> ST;
 
  626   RCP<Teuchos::FancyOStream> out = this->getOStream();
 
  627   Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
 
  629   if (!is_null(integrationControlStrategy_)) {
 
  630     integrationControlStrategy_->setOStream(out);
 
  631     integrationControlStrategy_->setVerbLevel(incrVerbLevel(verbLevel,-1));
 
  634   if (!is_null(integrationObserver_)) {
 
  635     integrationObserver_->setOStream(out);
 
  636     integrationObserver_->setVerbLevel(incrVerbLevel(verbLevel,-1));
 
  639   if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  640     *out << 
"\nEntering " << this->Describable::description()
 
  641          << 
"::advanceStepperToTime("<<advance_to_t<<
") ...\n";
 
  644   const int initCurrTimeStepIndex = currTimeStepIndex_;
 
  648   TimeRange<Scalar> currStepperTimeRange = stepper_->getTimeRange();
 
  651   bool return_val = 
true;
 
  653   while ( !currStepperTimeRange.isInRange(advance_to_t) ) {
 
  656     if (currTimeStepIndex_ >= maxNumTimeSteps_) {
 
  657       if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  660           << 
"\n*** NOTICE: currTimeStepIndex = "<<currTimeStepIndex_
 
  661           << 
" >= maxNumTimeSteps = "<<maxNumTimeSteps_
 
  662           << 
", halting time integration!" 
  668     if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) {
 
  669       *out << 
"\nTake step:  current_stepper_t = " 
  670            << currStepperTimeRange.upper()
 
  671            << 
", currTimeStepIndex = " << currTimeStepIndex_ << endl;
 
  680     if (stepCtrlInfoLast_.limitedByBreakPoint) {
 
  681       if ( stepCtrlInfoLast_.breakPointType == BREAK_POINT_TYPE_HARD ) {
 
  682         RYTHMOS_FUNC_TIME_MONITOR(
"Rythmos:DefaultIntegrator::restart");
 
  683         if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  684           *out << 
"\nAt a hard-breakpoint, restarting time integrator ...\n";
 
  688         if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  689           *out <<
"\nAt a soft-breakpoint, NOT restarting time integrator ...\n";
 
  699     bool foundAcceptableTimeStep = 
false;
 
  700     StepControlInfo<Scalar> stepCtrlInfo;
 
  705     while (!foundAcceptableTimeStep) {
 
  711       StepControlInfo<Scalar> trialStepCtrlInfo;
 
  713         RYTHMOS_FUNC_TIME_MONITOR(
 
  714           "Rythmos:DefaultIntegrator::advanceStepperToTime: getStepCtrl");
 
  715         if (!is_null(integrationControlStrategy_)) {
 
  719             integrationControlStrategy_->getNextStepControlInfo(
 
  720               *stepper_, stepCtrlInfoLast_, currTimeStepIndex_);
 
  724           trialStepCtrlInfo.stepType = STEP_TYPE_VARIABLE;
 
  725           trialStepCtrlInfo.stepSize = NL::max();
 
  730       if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  731         *out << 
"\nTrial step:\n";
 
  733         *out << trialStepCtrlInfo;
 
  737       if (trialStepCtrlInfo.stepSize < ST::zero()) {
 
  738         if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
 
  741             << 
"\n*** NOTICE: The IntegrationControlStrategy object " 
  742             << 
"return stepSize < 0.0, halting time integration!" 
  749       bool updatedTrialStepCtrlInfo = 
false;
 
  751         const Scalar finalTime = integrationTimeDomain_.upper();
 
  752         if (landOnFinalTime_ && trialStepCtrlInfo.stepSize
 
  753                                 + currStepperTimeRange.upper() > finalTime) {
 
  755           trialStepCtrlInfo.stepSize = finalTime - currStepperTimeRange.upper();
 
  756           updatedTrialStepCtrlInfo = 
true;
 
  758           if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  759             *out << 
"\nCutting trial step to "<< trialStepCtrlInfo.stepSize
 
  760                  << 
" to avoid stepping past final time ...\n";
 
  761           if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  763                  << 
"  integrationTimeDomain = " << integrationTimeDomain_<<
"\n" 
  764                  << 
"  currStepperTimeRange  = " << currStepperTimeRange <<
"\n";
 
  770       if ( updatedTrialStepCtrlInfo
 
  771         && includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
 
  773         *out << 
"\nUpdated trial step:\n";
 
  775         *out << trialStepCtrlInfo;
 
  783       if (!is_null(integrationObserver_))
 
  784         integrationObserver_->observeStartTimeStep(
 
  785             *stepper_, trialStepCtrlInfo, currTimeStepIndex_
 
  789       if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  790         if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE)
 
  791           *out << 
"\nTaking a variable time step with max step size = " 
  792                << trialStepCtrlInfo.stepSize << 
" ....\n";
 
  794           *out << 
"\nTaking a fixed time step of size = " 
  795                << trialStepCtrlInfo.stepSize << 
" ....\n";
 
  799       Scalar stepSizeTaken;
 
  801         RYTHMOS_FUNC_TIME_MONITOR(
 
  802           "Rythmos:Defaultintegrator::advancesteppertotime: takestep");
 
  803         stepSizeTaken = stepper_->takeStep(
 
  804           trialStepCtrlInfo.stepSize, trialStepCtrlInfo.stepType
 
  809       currStepperTimeRange = stepper_->getTimeRange();
 
  810       stepCtrlInfo = stepCtrlInfoTaken(trialStepCtrlInfo,stepSizeTaken);
 
  813       if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  814         *out << 
"\nStep actually taken:\n";
 
  816         *out << stepCtrlInfo;
 
  820       const bool timeStepFailed = (stepCtrlInfo.stepSize <= ST::zero());
 
  821       if (timeStepFailed && includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM)) {
 
  822         *out << 
"\nWARNING: timeStep = "<<trialStepCtrlInfo.stepSize
 
  827       if (timeStepFailed) {
 
  828         if (nonnull(integrationObserver_))
 
  829           integrationObserver_->observeFailedTimeStep(
 
  830             *stepper_, stepCtrlInfo, currTimeStepIndex_
 
  834         const bool handlesFailedTimeSteps =
 
  835           nonnull(integrationControlStrategy_) &&
 
  836           integrationControlStrategy_->handlesFailedTimeSteps();
 
  837         if (handlesFailedTimeSteps)
 
  840           if (integrationControlStrategy_->resetForFailedTimeStep(
 
  841                 *stepper_, stepCtrlInfoLast_,
 
  842                  currTimeStepIndex_, trialStepCtrlInfo) )
 
  844             if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  845               *out << 
"\nThe IntegrationControlStrategy object indicated that" 
  846                 << 
" it would like to suggest another timestep!\n";
 
  855             if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  856               *out << 
"\nThe IntegrationControlStrategy object could not " 
  857                    << 
"suggest a better time step!  Allowing to fail the " 
  866       if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE) {
 
  867         TEUCHOS_TEST_FOR_EXCEPTION(
 
  868           stepSizeTaken < ST::zero(), std::logic_error,
 
  869           "Error, stepper took negative step of dt = " << stepSizeTaken << 
"!\n" 
  871         TEUCHOS_TEST_FOR_EXCEPTION(
 
  872           stepSizeTaken > trialStepCtrlInfo.stepSize, std::logic_error,
 
  873           "Error, stepper took step of dt = " << stepSizeTaken
 
  874           << 
" > max step size of = " << trialStepCtrlInfo.stepSize << 
"!\n" 
  878         TEUCHOS_TEST_FOR_EXCEPTION(
 
  879           stepSizeTaken != trialStepCtrlInfo.stepSize, std::logic_error,
 
  880           "Error, stepper took step of dt = " << stepSizeTaken
 
  881           << 
" when asked to take step of dt = " << trialStepCtrlInfo.stepSize
 
  886       foundAcceptableTimeStep = 
true;
 
  889       if (!is_null(trailingInterpBuffer_)) {
 
  890         interpBufferAppender_->append(*stepper_,currStepperTimeRange,
 
  891           trailingInterpBuffer_.ptr() );
 
  900         RYTHMOS_FUNC_TIME_MONITOR(
 
  901           "Rythmos:DefaultIntegrator::advanceStepperToTime: output");
 
  904         if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
 
  905           StepStatus<Scalar> stepStatus = stepper_->getStepStatus();
 
  906           *out << 
"\nTime point reached = " << stepStatus.time << endl;
 
  907           *out << 
"\nstepStatus:\n" << stepStatus;
 
  908           if ( includesVerbLevel(verbLevel,Teuchos::VERB_EXTREME) ) {
 
  909             RCP<const Thyra::VectorBase<Scalar> >
 
  910               solution = stepStatus.solution,
 
  911               solutionDot = stepStatus.solutionDot;
 
  912             if (!is_null(solution))
 
  913               *out << 
"\nsolution = \n" 
  914                    << Teuchos::describe(*solution,verbLevel);
 
  915             if (!is_null(solutionDot))
 
  916               *out << 
"\nsolutionDot = \n" 
  917                    << Teuchos::describe(*solutionDot,verbLevel);
 
  922         if (!is_null(integrationObserver_))
 
  923           integrationObserver_->observeCompletedTimeStep(
 
  924             *stepper_, stepCtrlInfo, currTimeStepIndex_
 
  935     stepCtrlInfoLast_ = stepCtrlInfo;
 
  936     ++currTimeStepIndex_;
 
  940   if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
 
  941     *out << 
"\nNumber of steps taken in this call to " 
  942          << 
"advanceStepperToTime(...) = " 
  943          << (currTimeStepIndex_ - initCurrTimeStepIndex) << endl
 
  944          << 
"\nLeaving " << this->Describable::description()
 
  945          << 
"::advanceStepperToTime("<<advance_to_t<<
") ...\n";
 
  957 #define RYTHMOS_DEFAULT_INTEGRATOR_INSTANT(SCALAR) \ 
  959   template class DefaultIntegrator< SCALAR >; \ 
  961   template RCP<DefaultIntegrator< SCALAR > > \ 
  962   defaultIntegrator(); \ 
  964   template RCP<DefaultIntegrator< SCALAR > > \ 
  966     const RCP<IntegrationControlStrategyBase< SCALAR > > &integrationControlStrategy, \ 
  967     const RCP<IntegrationObserverBase< SCALAR > > &integrationObserver \ 
  970   template RCP<DefaultIntegrator< SCALAR > > \ 
  971   controlledDefaultIntegrator( \ 
  972     const RCP<IntegrationControlStrategyBase< SCALAR > > &integrationControlStrategy \ 
  975   template RCP<DefaultIntegrator< SCALAR > > \ 
  976   observedDefaultIntegrator( \ 
  977     const RCP<IntegrationObserverBase< SCALAR > > &integrationObserver \ 
  983 #endif //RYTHMOS_DEFAULT_INTEGRATOR_DEF_HPP 
RCP< const ParameterList > getValidParameters() const 
 
RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const 
 
void getPoints(const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec) const 
 
RCP< StepperBase< Scalar > > unSetStepper()
 
RCP< const InterpolationBufferAppenderBase< Scalar > > getInterpolationBufferAppender()
 
RCP< InterpolationBufferBase< Scalar > > unSetTrailingInterpolationBuffer()
 
Base class for defining stepper functionality. 
 
RCP< IntegratorBase< Scalar > > cloneIntegrator() const 
 
RCP< const StepperBase< Scalar > > getStepper() const 
 
RCP< IntegrationControlStrategyBase< Scalar > > getNonconstIntegrationControlStrategy()
 
Simple struct to aggregate integration/stepper control information. 
 
void getNodes(Array< Scalar > *time_vec) const 
 
Base class for strategy objects that control integration by selecting step sizes for a stepper...
 
RCP< const IntegrationControlStrategyBase< Scalar > > getIntegrationControlStrategy() const 
 
TimeRange< Scalar > getFwdTimeRange() const 
 
Base class for strategy objects that observe and time integration by observing the stepper object...
 
RCP< DefaultIntegrator< Scalar > > defaultIntegrator(const RCP< IntegrationControlStrategyBase< Scalar > > &integrationControlStrategy, const RCP< IntegrationObserverBase< Scalar > > &integrationObserver)
 
RCP< InterpolationBufferAppenderBase< Scalar > > unSetInterpolationBufferAppender()
 
RCP< DefaultIntegrator< Scalar > > defaultIntegrator()
 
void setIntegrationControlStrategy(const RCP< IntegrationControlStrategyBase< Scalar > > &integrationControlStrategy)
 
A concrete subclass for IntegratorBase that allows a good deal of customization. 
 
void addPoints(const Array< Scalar > &time_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &x_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &xdot_vec)
 
void setTrailingInterpolationBuffer(const RCP< InterpolationBufferBase< Scalar > > &trailingInterpBuffer)
 
RCP< StepperBase< Scalar > > getNonconstStepper() const 
 
Base class for an interpolation buffer. 
 
RCP< const InterpolationBufferBase< Scalar > > getTrailingInterpolationBuffer() const 
 
RCP< DefaultIntegrator< Scalar > > controlledDefaultIntegrator(const RCP< IntegrationControlStrategyBase< Scalar > > &integrationControlStrategy)
 
void getFwdPoints(const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec)
 
Base class for strategy objects that append data from one InterplationBufferBase object to another...
 
void removeNodes(Array< Scalar > &time_vec)
 
void setInterpolationBufferAppender(const RCP< InterpolationBufferAppenderBase< Scalar > > &interpBufferAppender)
 
RCP< InterpolationBufferBase< Scalar > > getNonconstTrailingInterpolationBuffer()
 
TimeRange< Scalar > getTimeRange() const 
 
RCP< DefaultIntegrator< Scalar > > observedDefaultIntegrator(const RCP< IntegrationObserverBase< Scalar > > &integrationObserver)
 
RCP< InterpolationBufferAppenderBase< Scalar > > getNonconstInterpolationBufferAppender()
 
void setIntegrationObserver(const RCP< IntegrationObserverBase< Scalar > > &integrationObserver)
 
void setParameterList(RCP< ParameterList > const ¶mList)
 
void setStepper(const RCP< StepperBase< Scalar > > &stepper, const Scalar &finalTime, const bool landOnFinalTime=true)