00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef RYTHMOS_DEFAULT_INTEGRATOR_HPP
00030 #define RYTHMOS_DEFAULT_INTEGRATOR_HPP
00031
00032
00033 #include "Rythmos_IntegratorBase.hpp"
00034 #include "Rythmos_InterpolationBufferHelpers.hpp"
00035 #include "Rythmos_IntegrationControlStrategyBase.hpp"
00036 #include "Rythmos_IntegrationObserverBase.hpp"
00037 #include "Rythmos_InterpolationBufferAppenderBase.hpp"
00038 #include "Rythmos_PointwiseInterpolationBufferAppender.hpp"
00039 #include "Rythmos_StepperHelpers.hpp"
00040 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
00041 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
00042 #include "Teuchos_Assert.hpp"
00043 #include "Teuchos_as.hpp"
00044
00045
00046 namespace Rythmos {
00047
00048
00052 template<class Scalar>
00053 class DefaultIntegrator
00054 : virtual public IntegratorBase<Scalar>,
00055 virtual public Teuchos::ParameterListAcceptorDefaultBase
00056 {
00057 public:
00058
00060 typedef typename ScalarTraits<Scalar>::magnitudeType ScalarMag;
00061
00064
00066 DefaultIntegrator();
00067
00069 void setIntegrationControlStrategy(
00070 const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy
00071 );
00072
00074 void setIntegrationObserver(
00075 const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00076 );
00077
00079 void setInterpolationBufferAppender(
00080 const RCP<InterpolationBufferAppenderBase<Scalar> > &interpBufferAppender
00081 );
00082
00084 RCP<const InterpolationBufferAppenderBase<Scalar> >
00085 getInterpolationBufferAppender();
00086
00088 RCP<InterpolationBufferAppenderBase<Scalar> >
00089 unSetInterpolationBufferAppender();
00090
00092
00095
00097 void setParameterList(RCP<ParameterList> const& paramList);
00098
00100 RCP<const ParameterList> getValidParameters() const;
00101
00103
00106
00108 RCP<IntegratorBase<Scalar> > cloneIntegrator() const;
00109
00111 void setStepper(
00112 const RCP<StepperBase<Scalar> > &stepper,
00113 const Scalar &finalTime,
00114 const bool landOnFinalTime
00115 );
00116
00118 RCP<StepperBase<Scalar> > unSetStepper();
00119
00121 RCP<const StepperBase<Scalar> > getStepper() const;
00122
00124 bool acceptsTrailingInterpolationBuffer() const;
00125
00127 void setTrailingInterpolationBuffer(
00128 const RCP<InterpolationBufferBase<Scalar> > &trailingInterpBuffer
00129 );
00130
00132 RCP<InterpolationBufferBase<Scalar> >
00133 getTrailingInterpolationBuffer();
00134
00136 RCP<const InterpolationBufferBase<Scalar> >
00137 getTrailingInterpolationBuffer() const;
00138
00140 void getFwdPoints(
00141 const Array<Scalar>& time_vec,
00142 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00143 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00144 Array<ScalarMag>* accuracy_vec
00145 );
00146
00148 TimeRange<Scalar> getFwdTimeRange() const;
00149
00151
00154
00156 RCP<const Thyra::VectorSpaceBase<Scalar> > get_x_space() const;
00157
00159 void addPoints(
00160 const Array<Scalar>& time_vec,
00161 const Array<RCP<const Thyra::VectorBase<Scalar> > >& x_vec,
00162 const Array<RCP<const Thyra::VectorBase<Scalar> > >& xdot_vec
00163 );
00164
00166 void getPoints(
00167 const Array<Scalar>& time_vec,
00168 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00169 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00170 Array<ScalarMag>* accuracy_vec
00171 ) const;
00172
00174 TimeRange<Scalar> getTimeRange() const;
00175
00177 void getNodes(Array<Scalar>* time_vec) const;
00178
00180 void removeNodes(Array<Scalar>& time_vec);
00181
00183 int getOrder() const;
00184
00186
00187 private:
00188
00189
00190
00191
00192 RCP<IntegrationControlStrategyBase<Scalar> > integrationControlStrategy_;
00193 RCP<IntegrationObserverBase<Scalar> > integrationObserver_;
00194
00195 RCP<InterpolationBufferBase<Scalar> > trailingInterpBuffer_;
00196 RCP<InterpolationBufferAppenderBase<Scalar> > interpBufferAppender_;
00197
00198 RCP<StepperBase<Scalar> > stepper_;
00199 TimeRange<Scalar> integrationTimeDomain_;
00200 bool landOnFinalTime_;
00201
00202 int maxNumTimeSteps_;
00203
00204 int currTimeStepIndex_;
00205 StepControlInfo<Scalar> stepCtrlInfoLast_;
00206
00207 static const std::string maxNumTimeSteps_name_;
00208 static const int maxNumTimeSteps_default_;
00209
00210
00211
00212
00213 void finalizeSetup();
00214
00215 bool advanceStepperToTime( const Scalar& t );
00216
00217 };
00218
00219
00224 template<class Scalar>
00225 RCP<DefaultIntegrator<Scalar> >
00226 defaultIntegrator()
00227 {
00228 RCP<DefaultIntegrator<Scalar> >
00229 integrator = Teuchos::rcp(new DefaultIntegrator<Scalar>());
00230 return integrator;
00231 }
00232
00233
00238 template<class Scalar>
00239 RCP<DefaultIntegrator<Scalar> >
00240 defaultIntegrator(
00241 const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy,
00242 const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00243 )
00244 {
00245 RCP<DefaultIntegrator<Scalar> >
00246 integrator = Teuchos::rcp(new DefaultIntegrator<Scalar>());
00247 integrator->setIntegrationControlStrategy(integrationControlStrategy);
00248 integrator->setIntegrationObserver(integrationObserver);
00249 return integrator;
00250 }
00251
00252
00257 template<class Scalar>
00258 RCP<DefaultIntegrator<Scalar> >
00259 controlledDefaultIntegrator(
00260 const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy
00261 )
00262 {
00263 RCP<DefaultIntegrator<Scalar> >
00264 integrator = Teuchos::rcp(new DefaultIntegrator<Scalar>());
00265 integrator->setIntegrationControlStrategy(integrationControlStrategy);
00266 return integrator;
00267 }
00268
00269
00274 template<class Scalar>
00275 RCP<DefaultIntegrator<Scalar> >
00276 observedDefaultIntegrator(
00277 const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00278 )
00279 {
00280 RCP<DefaultIntegrator<Scalar> >
00281 integrator = Teuchos::rcp(new DefaultIntegrator<Scalar>());
00282 integrator->setIntegrationObserver(integrationObserver);
00283 return integrator;
00284 }
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300 template<class Scalar>
00301 const std::string
00302 DefaultIntegrator<Scalar>::maxNumTimeSteps_name_ = "Max Number Time Steps";
00303
00304 template<class Scalar>
00305 const int
00306 DefaultIntegrator<Scalar>::maxNumTimeSteps_default_ = 10000;
00307
00308
00309
00310
00311
00312 template<class Scalar>
00313 DefaultIntegrator<Scalar>::DefaultIntegrator()
00314 :landOnFinalTime_(true),
00315 maxNumTimeSteps_(maxNumTimeSteps_default_),
00316 currTimeStepIndex_(-1)
00317 {}
00318
00319
00320 template<class Scalar>
00321 void DefaultIntegrator<Scalar>::setIntegrationControlStrategy(
00322 const RCP<IntegrationControlStrategyBase<Scalar> > &integrationControlStrategy
00323 )
00324 {
00325 #ifdef TEUCHOS_DEBUG
00326 TEST_FOR_EXCEPT(is_null(integrationControlStrategy));
00327 #endif
00328 integrationControlStrategy_ = integrationControlStrategy;
00329 }
00330
00331
00332 template<class Scalar>
00333 void DefaultIntegrator<Scalar>::setIntegrationObserver(
00334 const RCP<IntegrationObserverBase<Scalar> > &integrationObserver
00335 )
00336 {
00337 #ifdef TEUCHOS_DEBUG
00338 TEST_FOR_EXCEPT(is_null(integrationObserver));
00339 #endif
00340 integrationObserver_ = integrationObserver;
00341 }
00342
00343
00344 template<class Scalar>
00345 void DefaultIntegrator<Scalar>::setInterpolationBufferAppender(
00346 const RCP<InterpolationBufferAppenderBase<Scalar> > &interpBufferAppender
00347 )
00348 {
00349 interpBufferAppender_ = interpBufferAppender.assert_not_null();
00350 }
00351
00352
00353 template<class Scalar>
00354 RCP<const InterpolationBufferAppenderBase<Scalar> >
00355 DefaultIntegrator<Scalar>::getInterpolationBufferAppender()
00356 {
00357 return interpBufferAppender_;
00358 }
00359
00360
00361 template<class Scalar>
00362 RCP<InterpolationBufferAppenderBase<Scalar> >
00363 DefaultIntegrator<Scalar>::unSetInterpolationBufferAppender()
00364 {
00365 RCP<InterpolationBufferAppenderBase<Scalar> > InterpBufferAppender;
00366 std::swap( InterpBufferAppender, interpBufferAppender_ );
00367 return InterpBufferAppender;
00368 }
00369
00370
00371
00372
00373
00374 template<class Scalar>
00375 void DefaultIntegrator<Scalar>::setParameterList(
00376 RCP<ParameterList> const& paramList
00377 )
00378 {
00379 TEST_FOR_EXCEPT(is_null(paramList));
00380 paramList->validateParameters(*getValidParameters());
00381 this->setMyParamList(paramList);
00382 maxNumTimeSteps_ = paramList->get(
00383 maxNumTimeSteps_name_, maxNumTimeSteps_default_);
00384 Teuchos::readVerboseObjectSublist(&*paramList,this);
00385 }
00386
00387
00388 template<class Scalar>
00389 RCP<const ParameterList>
00390 DefaultIntegrator<Scalar>::getValidParameters() const
00391 {
00392 static RCP<const ParameterList> validPL;
00393 if (is_null(validPL) ) {
00394 RCP<ParameterList> pl = Teuchos::parameterList();
00395 pl->set(maxNumTimeSteps_name_, maxNumTimeSteps_default_,
00396 "Set the maximum number of integration time-steps allowed."
00397 );
00398 Teuchos::setupVerboseObjectSublist(&*pl);
00399 validPL = pl;
00400 }
00401 return validPL;
00402 }
00403
00404
00405
00406
00407
00408 template<class Scalar>
00409 RCP<IntegratorBase<Scalar> >
00410 DefaultIntegrator<Scalar>::cloneIntegrator() const
00411 {
00412 RCP<DefaultIntegrator<Scalar> >
00413 newIntegrator = Teuchos::rcp(new DefaultIntegrator<Scalar>());
00414
00415 newIntegrator->stepper_ = Teuchos::null;
00416 const RCP<const ParameterList> paramList = this->getParameterList();
00417 if (!is_null(paramList)) {
00418 newIntegrator->setParameterList(Teuchos::parameterList(*paramList));
00419 }
00420 if (!is_null(integrationControlStrategy_)) {
00421 newIntegrator->integrationControlStrategy_ =
00422 integrationControlStrategy_->cloneIntegrationControlStrategy().assert_not_null();
00423 }
00424 if (!is_null(integrationObserver_)) {
00425 newIntegrator->integrationObserver_ =
00426 integrationObserver_->cloneIntegrationObserver().assert_not_null();
00427 }
00428 if (!is_null(trailingInterpBuffer_)) {
00429 TEST_FOR_EXCEPT_MSG(true,"ToDo: Clone the trailing interpolation buffer");
00430 }
00431 if (!is_null(interpBufferAppender_)) {
00432 TEST_FOR_EXCEPT_MSG(true,"ToDo: Clone the trailing interpolation buffer appender");
00433 }
00434 return newIntegrator;
00435 }
00436
00437
00438 template<class Scalar>
00439 void DefaultIntegrator<Scalar>::setStepper(
00440 const RCP<StepperBase<Scalar> > &stepper,
00441 const Scalar &finalTime,
00442 const bool landOnFinalTime
00443 )
00444 {
00445 typedef Teuchos::ScalarTraits<Scalar> ST;
00446 TEST_FOR_EXCEPT(is_null(stepper));
00447 TEST_FOR_EXCEPT( finalTime <= stepper->getTimeRange().lower() );
00448 TEUCHOS_ASSERT( stepper->getTimeRange().length() == ST::zero() );
00449
00450 stepper_ = stepper;
00451 integrationTimeDomain_ = timeRange(stepper_->getTimeRange().lower(), finalTime);
00452 landOnFinalTime_ = landOnFinalTime;
00453 currTimeStepIndex_ = 0;
00454 stepCtrlInfoLast_ = StepControlInfo<Scalar>();
00455 if (!is_null(integrationControlStrategy_))
00456 integrationControlStrategy_->resetIntegrationControlStrategy(
00457 integrationTimeDomain_
00458 );
00459 if (!is_null(integrationObserver_))
00460 integrationObserver_->resetIntegrationObserver(
00461 integrationTimeDomain_
00462 );
00463 }
00464
00465
00466 template<class Scalar>
00467 RCP<StepperBase<Scalar> > DefaultIntegrator<Scalar>::unSetStepper()
00468 {
00469 RCP<StepperBase<Scalar> > stepper_temp = stepper_;
00470 stepper_ = Teuchos::null;
00471 return(stepper_temp);
00472 }
00473
00474
00475 template<class Scalar>
00476 RCP<const StepperBase<Scalar> > DefaultIntegrator<Scalar>::getStepper() const
00477 {
00478 return(stepper_);
00479 }
00480
00481
00482 template<class Scalar>
00483 bool DefaultIntegrator<Scalar>::acceptsTrailingInterpolationBuffer() const
00484 {
00485 return true;
00486 }
00487
00488
00489 template<class Scalar>
00490 void DefaultIntegrator<Scalar>::setTrailingInterpolationBuffer(
00491 const RCP<InterpolationBufferBase<Scalar> > &trailingInterpBuffer
00492 )
00493 {
00494 trailingInterpBuffer_ = trailingInterpBuffer.assert_not_null();
00495 }
00496
00497
00498 template<class Scalar>
00499 RCP<InterpolationBufferBase<Scalar> >
00500 DefaultIntegrator<Scalar>::getTrailingInterpolationBuffer()
00501 {
00502 return trailingInterpBuffer_;
00503 }
00504
00505
00506 template<class Scalar>
00507 RCP<const InterpolationBufferBase<Scalar> >
00508 DefaultIntegrator<Scalar>::getTrailingInterpolationBuffer() const
00509 {
00510 return trailingInterpBuffer_;
00511 }
00512
00513
00514 template<class Scalar>
00515 void DefaultIntegrator<Scalar>::getFwdPoints(
00516 const Array<Scalar>& time_vec,
00517 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00518 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00519 Array<ScalarMag>* accuracy_vec
00520 )
00521 {
00522
00523 #ifdef ENABLE_RYTHMOS_TIMERS
00524 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints");
00525 #endif
00526
00527 using Teuchos::incrVerbLevel;
00528 using Teuchos::Describable;
00529 typedef Teuchos::ScalarTraits<Scalar> ST;
00530 typedef InterpolationBufferBase<Scalar> IBB;
00531 typedef Teuchos::VerboseObjectTempState<IBB> VOTSIBB;
00532
00533 finalizeSetup();
00534
00535 RCP<Teuchos::FancyOStream> out = this->getOStream();
00536 Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
00537 Teuchos::OSTab tab(out);
00538 VOTSIBB stepper_outputTempState(stepper_,out,incrVerbLevel(verbLevel,-1));
00539
00540 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00541 *out << "\nEntering " << this->Describable::description() << "::getFwdPoints(...) ...\n"
00542 << "\nStepper: " << Teuchos::describe(*stepper_,verbLevel);
00543
00544 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
00545 *out << "\nRequested time points: " << Teuchos::toString(time_vec) << "\n";
00546
00547
00548
00549
00550
00551 const int numTimePoints = time_vec.size();
00552
00553
00554 assertTimePointsAreSorted(time_vec);
00555 TEST_FOR_EXCEPT(accuracy_vec!=0);
00556
00557
00558 if (x_vec)
00559 x_vec->resize(numTimePoints);
00560 if (xdot_vec)
00561 xdot_vec->resize(numTimePoints);
00562
00563
00564
00565
00566 int nextTimePointIndex = 0;
00567
00568 assertNoTimePointsBeforeCurrentTimeRange(*this,time_vec,nextTimePointIndex);
00569
00570
00571
00572
00573
00574 {
00575 #ifdef ENABLE_RYTHMOS_TIMERS
00576 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: getPoints");
00577 #endif
00578
00579 getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex);
00580 }
00581
00582
00583
00584
00585
00586
00587 while ( nextTimePointIndex < numTimePoints ) {
00588
00589
00590
00591 const Scalar t = time_vec[nextTimePointIndex];
00592 bool advanceStepperToTimeSucceeded = false;
00593 {
00594 #ifdef ENABLE_RYTHMOS_TIMERS
00595 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: advanceStepperToTime");
00596 #endif
00597 advanceStepperToTimeSucceeded= advanceStepperToTime(t);
00598 }
00599 TEST_FOR_EXCEPTION(
00600 !advanceStepperToTimeSucceeded, Exceptions::GetFwdPointsFailed,
00601 this->description() << "\n\n"
00602 "Error: The integration failed to get to time " << t << " and only achieved\n"
00603 "getting to " << stepper_->getTimeRange().upper() << "!"
00604 );
00605
00606
00607 {
00608 #ifdef ENABLE_RYTHMOS_TIMERS
00609 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: getPoints (fwd)");
00610 #endif
00611 getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex);
00612 }
00613
00614 }
00615
00616 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00617 *out << "\nLeaving " << this->Describable::description() << "::getFwdPoints(...) ...\n";
00618
00619 }
00620
00621
00622 template<class Scalar>
00623 TimeRange<Scalar>
00624 DefaultIntegrator<Scalar>::getFwdTimeRange() const
00625 {
00626 return timeRange(
00627 stepper_->getTimeRange().lower(),
00628 integrationTimeDomain_.upper()
00629 );
00630 }
00631
00632
00633
00634
00635
00636 template<class Scalar>
00637 RCP<const Thyra::VectorSpaceBase<Scalar> >
00638 DefaultIntegrator<Scalar>::get_x_space() const
00639 {
00640 return stepper_->get_x_space();
00641 }
00642
00643
00644 template<class Scalar>
00645 void DefaultIntegrator<Scalar>::addPoints(
00646 const Array<Scalar>& time_vec,
00647 const Array<RCP<const Thyra::VectorBase<Scalar> > >& x_vec,
00648 const Array<RCP<const Thyra::VectorBase<Scalar> > >& xdot_vec
00649 )
00650 {
00651 stepper_->addPoints(time_vec,x_vec,xdot_vec);
00652 }
00653
00654
00655 template<class Scalar>
00656 void DefaultIntegrator<Scalar>::getPoints(
00657 const Array<Scalar>& time_vec,
00658 Array<RCP<const Thyra::VectorBase<Scalar> > >* x_vec,
00659 Array<RCP<const Thyra::VectorBase<Scalar> > >* xdot_vec,
00660 Array<ScalarMag>* accuracy_vec
00661 ) const
00662 {
00663
00664 stepper_->getPoints(time_vec,x_vec,xdot_vec,accuracy_vec);
00665 }
00666
00667
00668 template<class Scalar>
00669 TimeRange<Scalar> DefaultIntegrator<Scalar>::getTimeRange() const
00670 {
00671 return stepper_->getTimeRange();
00672 }
00673
00674
00675 template<class Scalar>
00676 void DefaultIntegrator<Scalar>::getNodes(Array<Scalar>* time_vec) const
00677 {
00678 stepper_->getNodes(time_vec);
00679 }
00680
00681
00682 template<class Scalar>
00683 void DefaultIntegrator<Scalar>::removeNodes(Array<Scalar>& time_vec)
00684 {
00685 stepper_->removeNodes(time_vec);
00686 }
00687
00688
00689 template<class Scalar>
00690 int DefaultIntegrator<Scalar>::getOrder() const
00691 {
00692 return stepper_->getOrder();
00693 }
00694
00695
00696
00697
00698
00699 template<class Scalar>
00700 void DefaultIntegrator<Scalar>::finalizeSetup()
00701 {
00702 if (!is_null(trailingInterpBuffer_) && is_null(interpBufferAppender_))
00703 interpBufferAppender_ = pointwiseInterpolationBufferAppender<Scalar>();
00704
00705 }
00706
00707
00708 template<class Scalar>
00709 bool DefaultIntegrator<Scalar>::advanceStepperToTime( const Scalar& advance_to_t )
00710 {
00711
00712 #ifdef ENABLE_RYTHMOS_TIMERS
00713 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime");
00714 #endif
00715
00716 using std::endl;
00717 typedef std::numeric_limits<Scalar> NL;
00718 using Teuchos::incrVerbLevel;
00719 using Teuchos::Describable;
00720 using Teuchos::OSTab;
00721 typedef Teuchos::ScalarTraits<Scalar> ST;
00722
00723 RCP<Teuchos::FancyOStream> out = this->getOStream();
00724 Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
00725 Teuchos::OSTab tab(out);
00726
00727 if (!is_null(integrationControlStrategy_)) {
00728 integrationControlStrategy_->setOStream(out);
00729 integrationControlStrategy_->setVerbLevel(incrVerbLevel(verbLevel,-1));
00730 }
00731
00732 if (!is_null(integrationObserver_)) {
00733 integrationObserver_->setOStream(out);
00734 integrationObserver_->setVerbLevel(incrVerbLevel(verbLevel,-1));
00735 }
00736
00737 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00738 *out << "\nEntering " << this->Describable::description()
00739 << "::advanceStepperToTime("<<advance_to_t<<") ...\n";
00740
00741
00742 const int initCurrTimeStepIndex = currTimeStepIndex_;
00743
00744
00745
00746 TimeRange<Scalar> currStepperTimeRange = stepper_->getTimeRange();
00747
00748
00749 bool return_val = true;
00750
00751 while ( !currStepperTimeRange.isInRange(advance_to_t) ) {
00752
00753
00754 if (currTimeStepIndex_ >= maxNumTimeSteps_) {
00755 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00756 *out
00757 << "\n***"
00758 << "\n*** NOTICE: currTimeStepIndex = "<<currTimeStepIndex_
00759 << " >= maxNumTimeSteps = "<<maxNumTimeSteps_<< ", halting time integration!"
00760 << "\n***\n";
00761 return_val = false;
00762 break;
00763 }
00764
00765 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00766 *out << "\nTake step: current_stepper_t = " << currStepperTimeRange.upper()
00767 << ", currTimeStepIndex = " << currTimeStepIndex_ << endl;
00768 Teuchos::OSTab tab(out);
00769
00770
00771
00772
00773
00774 if (stepCtrlInfoLast_.limitedByBreakPoint) {
00775 if ( stepCtrlInfoLast_.breakPointType == BREAK_POINT_TYPE_HARD ) {
00776 #ifdef ENABLE_RYTHMOS_TIMERS
00777 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::restart");
00778 #endif
00779 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00780 *out << "\nAt a hard-breakpoint, restarting time integrator ...\n";
00781 restart(&*stepper_);
00782 }
00783 else {
00784 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00785 *out << "\nAt a soft-breakpoint, NOT restarting time integrator ...\n";
00786 }
00787 }
00788
00789
00790
00791
00792
00793 StepControlInfo<Scalar> trialStepCtrlInfo;
00794 {
00795 #ifdef ENABLE_RYTHMOS_TIMERS
00796 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: getStepCtrl");
00797 #endif
00798 if (!is_null(integrationControlStrategy_)) {
00799
00800
00801 trialStepCtrlInfo = integrationControlStrategy_->getNextStepControlInfo(
00802 *stepper_, stepCtrlInfoLast_, currTimeStepIndex_
00803 );
00804 }
00805 else {
00806
00807 trialStepCtrlInfo.stepType = STEP_TYPE_VARIABLE;
00808 trialStepCtrlInfo.stepSize = NL::max();
00809 }
00810 }
00811
00812
00813 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
00814 *out << "\nTrial step:\n";
00815 OSTab tab(out);
00816 *out << trialStepCtrlInfo;
00817 }
00818
00819
00820 if (trialStepCtrlInfo.stepSize < ST::zero()) {
00821 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
00822 *out
00823 << "\n***"
00824 << "\n*** NOTICE: The IntegrationControlStrategy object return stepSize < 0.0, halting time integration!"
00825 << "\n***\n";
00826 return_val = false;
00827 break;
00828 }
00829
00830
00831 bool updatedTrialStepCtrlInfo = false;
00832 {
00833 const Scalar finalTime = integrationTimeDomain_.upper();
00834 if (landOnFinalTime_ && trialStepCtrlInfo.stepSize + currStepperTimeRange.upper() > finalTime) {
00835 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00836 *out << "\nCutting trial step to avoid stepping past final time ...\n";
00837 trialStepCtrlInfo.stepSize = finalTime - currStepperTimeRange.upper();
00838 updatedTrialStepCtrlInfo = true;
00839 }
00840 }
00841
00842
00843 if ( updatedTrialStepCtrlInfo
00844 && includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) )
00845 {
00846 *out << "\nUpdated trial step:\n";
00847 OSTab tab(out);
00848 *out << trialStepCtrlInfo;
00849 }
00850
00851
00852
00853
00854
00855
00856 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
00857 if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE)
00858 *out << "\nTaking a variable time step with max step size = "
00859 << trialStepCtrlInfo.stepSize << " ....\n";
00860 else
00861 *out << "\nTaking a fixed time step of size = "
00862 << trialStepCtrlInfo.stepSize << " ....\n";
00863 }
00864
00865
00866 Scalar stepSizeTaken;
00867 {
00868 #ifdef ENABLE_RYTHMOS_TIMERS
00869 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: takeStep");
00870 #endif
00871 stepSizeTaken = stepper_->takeStep(
00872 trialStepCtrlInfo.stepSize, trialStepCtrlInfo.stepType
00873 );
00874 }
00875
00876
00877 if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE) {
00878 TEST_FOR_EXCEPTION(
00879 stepSizeTaken < ST::zero(), std::logic_error,
00880 "Error, stepper took negative step of dt = " << stepSizeTaken << "!\n"
00881 );
00882 TEST_FOR_EXCEPTION(
00883 stepSizeTaken > trialStepCtrlInfo.stepSize, std::logic_error,
00884 "Error, stepper took step of dt = " << stepSizeTaken
00885 << " > max step size of = " << trialStepCtrlInfo.stepSize << "!\n"
00886 );
00887 }
00888 else {
00889 TEST_FOR_EXCEPTION(
00890 stepSizeTaken != trialStepCtrlInfo.stepSize, std::logic_error,
00891 "Error, stepper took step of dt = " << stepSizeTaken
00892 << " when asked to take step of dt = " << trialStepCtrlInfo.stepSize << "\n"
00893 );
00894 }
00895
00896
00897 currStepperTimeRange = stepper_->getTimeRange();
00898 const StepControlInfo<Scalar> stepCtrlInfo =
00899 stepCtrlInfoTaken(trialStepCtrlInfo,stepSizeTaken);
00900
00901
00902 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
00903 *out << "\nStep actually taken:\n";
00904 OSTab tab(out);
00905 *out << stepCtrlInfo;
00906 }
00907
00908
00909 if (!is_null(trailingInterpBuffer_)) {
00910 interpBufferAppender_->append(*stepper_,currStepperTimeRange,
00911 trailingInterpBuffer_.ptr() );
00912 }
00913
00914
00915
00916
00917
00918 {
00919
00920 #ifdef ENABLE_RYTHMOS_TIMERS
00921 TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: output");
00922 #endif
00923
00924
00925 if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) {
00926 StepStatus<Scalar> stepStatus = stepper_->getStepStatus();
00927 *out << "\nTime point reached = " << stepStatus.time << endl;
00928 *out << "\nstepStatus:\n" << stepStatus;
00929 if ( includesVerbLevel(verbLevel,Teuchos::VERB_EXTREME) ) {
00930 RCP<const Thyra::VectorBase<Scalar> >
00931 solution = stepStatus.solution,
00932 solutionDot = stepStatus.solutionDot;
00933 if (!is_null(solution))
00934 *out << "\nsolution = \n" << Teuchos::describe(*solution,verbLevel);
00935 if (!is_null(solutionDot))
00936 *out << "\nsolutionDot = \n" << Teuchos::describe(*solutionDot,verbLevel);
00937 }
00938 }
00939
00940
00941 if (!is_null(integrationObserver_))
00942 integrationObserver_->observeCompletedTimeStep(
00943 *stepper_, stepCtrlInfo, currTimeStepIndex_
00944 );
00945
00946 }
00947
00948
00949
00950
00951
00952 stepCtrlInfoLast_ = stepCtrlInfo;
00953 ++currTimeStepIndex_;
00954
00955 }
00956
00957 if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) )
00958 *out << "\nNumber of steps taken in this call to advanceStepperToTime(...) = "
00959 << (currTimeStepIndex_ - initCurrTimeStepIndex) << endl
00960 << "\nLeaving" << this->Describable::description()
00961 << "::advanceStepperToTime("<<advance_to_t<<") ...\n";
00962
00963 return return_val;
00964
00965 }
00966
00967
00968 }
00969
00970
00971 #endif //RYTHMOS_DEFAULT_INTEGRATOR_HPP