bool DefaultIntegrator<Scalar>::advanceStepperToTime( const Scalar& advance_to_t ) { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR_DIFF("Rythmos:DefaultIntegrator::advanceStepperToTime", TopLevel); #endif using std::endl; typedef std::numeric_limits<Scalar> NL; using Teuchos::incrVerbLevel; #ifndef _MSC_VER using Teuchos::Describable; #endif using Teuchos::OSTab; typedef Teuchos::ScalarTraits<Scalar> ST; RCP<Teuchos::FancyOStream> out = this->getOStream(); Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel(); if (!is_null(integrationControlStrategy_)) { integrationControlStrategy_->setOStream(out); integrationControlStrategy_->setVerbLevel(incrVerbLevel(verbLevel,-1)); } if (!is_null(integrationObserver_)) { integrationObserver_->setOStream(out); integrationObserver_->setVerbLevel(incrVerbLevel(verbLevel,-1)); } if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nEntering " << this->Describable::description() << "::advanceStepperToTime("<<advance_to_t<<") ...\n"; // Remember what timestep index we are on so we can report it later const int initCurrTimeStepIndex = currTimeStepIndex_; // Take steps until we the requested time is reached (or passed) TimeRange<Scalar> currStepperTimeRange = stepper_->getTimeRange(); // Start by assume we can reach the time advance_to_t bool return_val = true; while ( !currStepperTimeRange.isInRange(advance_to_t) ) { // Halt immediatly if exceeded max iterations if (currTimeStepIndex_ >= maxNumTimeSteps_) { if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\n***" << "\n*** NOTICE: currTimeStepIndex = "<<currTimeStepIndex_ << " >= maxNumTimeSteps = "<<maxNumTimeSteps_<< ", halting time integration!" << "\n***\n"; return_val = false; break; // Exit the loop immediately! } if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nTake step: current_stepper_t = " << currStepperTimeRange.upper() << ", currTimeStepIndex = " << currTimeStepIndex_ << endl; // // A) Reinitialize if a hard breakpoint was reached on the last time step // if (stepCtrlInfoLast_.limitedByBreakPoint) { if ( stepCtrlInfoLast_.breakPointType == BREAK_POINT_TYPE_HARD ) { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::restart"); #endif if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nAt a hard-breakpoint, restarting time integrator ...\n"; restart(&*stepper_); } else { if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nAt a soft-breakpoint, NOT restarting time integrator ...\n"; } } // // B) Get the trial step control info // StepControlInfo<Scalar> trialStepCtrlInfo; { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: getStepCtrl"); #endif if (!is_null(integrationControlStrategy_)) { // Let an external strategy object determine the step size and type. // Note that any breakpoint info is also related through this call. trialStepCtrlInfo = integrationControlStrategy_->getNextStepControlInfo( *stepper_, stepCtrlInfoLast_, currTimeStepIndex_ ); } else { // Take a variable step if we have no control strategy trialStepCtrlInfo.stepType = STEP_TYPE_VARIABLE; trialStepCtrlInfo.stepSize = NL::max(); } } // Print the initial trial step if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) { *out << "\nTrial step:\n"; OSTab tab(out); *out << trialStepCtrlInfo; } // Halt immediately if we where told to do so if (trialStepCtrlInfo.stepSize < ST::zero()) { if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) *out << "\n***" << "\n*** NOTICE: The IntegrationControlStrategy object return stepSize < 0.0, halting time integration!" << "\n***\n"; return_val = false; break; // Exit the loop immediately! } // Make sure we don't step past the final time if asked not to bool updatedTrialStepCtrlInfo = false; { const Scalar finalTime = integrationTimeDomain_.upper(); if (landOnFinalTime_ && trialStepCtrlInfo.stepSize + currStepperTimeRange.upper() > finalTime) { if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nCutting trial step to avoid stepping past final time ...\n"; trialStepCtrlInfo.stepSize = finalTime - currStepperTimeRange.upper(); updatedTrialStepCtrlInfo = true; } } // Print the modified trial step if ( updatedTrialStepCtrlInfo && includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) { *out << "\nUpdated trial step:\n"; OSTab tab(out); *out << trialStepCtrlInfo; } // // C) Take the step // // Print step type and size if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) { if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE) *out << "\nTaking a variable time step with max step size = " << trialStepCtrlInfo.stepSize << " ....\n"; else *out << "\nTaking a fixed time step of size = " << trialStepCtrlInfo.stepSize << " ....\n"; } // Take step Scalar stepSizeTaken; { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: takeStep"); #endif stepSizeTaken = stepper_->takeStep( trialStepCtrlInfo.stepSize, trialStepCtrlInfo.stepType ); } // Validate step taken if (trialStepCtrlInfo.stepType == STEP_TYPE_VARIABLE) { TEST_FOR_EXCEPTION( stepSizeTaken < ST::zero(), std::logic_error, "Error, stepper took negative step of dt = " << stepSizeTaken << "!\n" ); TEST_FOR_EXCEPTION( stepSizeTaken > trialStepCtrlInfo.stepSize, std::logic_error, "Error, stepper took step of dt = " << stepSizeTaken << " > max step size of = " << trialStepCtrlInfo.stepSize << "!\n" ); } else { // STEP_TYPE_FIXED TEST_FOR_EXCEPTION( stepSizeTaken != trialStepCtrlInfo.stepSize, std::logic_error, "Error, stepper took step of dt = " << stepSizeTaken << " when asked to take step of dt = " << trialStepCtrlInfo.stepSize << "\n" ); } // Update info about this step currStepperTimeRange = stepper_->getTimeRange(); const StepControlInfo<Scalar> stepCtrlInfo = stepCtrlInfoTaken(trialStepCtrlInfo,stepSizeTaken); // Print the step actually taken if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) { *out << "\nStep actually taken:\n"; OSTab tab(out); *out << stepCtrlInfo; } // Append the trailing interpolation buffer (if defined) if (!is_null(trailingInterpBuffer_)) { interpBufferAppender_->append(*stepper_,currStepperTimeRange, trailingInterpBuffer_.ptr() ); } // // D) Output info about step // { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::advanceStepperToTime: output"); #endif // Print our own brief output if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) { StepStatus<Scalar> stepStatus = stepper_->getStepStatus(); *out << "\nTime point reached = " << stepStatus.time << endl; *out << "\nstepStatus:\n" << stepStatus; if ( includesVerbLevel(verbLevel,Teuchos::VERB_EXTREME) ) { RCP<const Thyra::VectorBase<Scalar> > solution = stepStatus.solution, solutionDot = stepStatus.solutionDot; if (!is_null(solution)) *out << "\nsolution = \n" << Teuchos::describe(*solution,verbLevel); if (!is_null(solutionDot)) *out << "\nsolutionDot = \n" << Teuchos::describe(*solutionDot,verbLevel); } } // Output to the observer if (!is_null(integrationObserver_)) integrationObserver_->observeCompletedTimeStep( *stepper_, stepCtrlInfo, currTimeStepIndex_ ); } // // E) Update info for next time step // stepCtrlInfoLast_ = stepCtrlInfo; ++currTimeStepIndex_; } if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nNumber of steps taken in this call to advanceStepperToTime(...) = " << (currTimeStepIndex_ - initCurrTimeStepIndex) << endl << "\nLeaving" << this->Describable::description() << "::advanceStepperToTime("<<advance_to_t<<") ...\n"; return return_val; }
Scalar ImplicitRKStepper<Scalar>::takeStep(Scalar dt, StepSizeType stepSizeType) { using Teuchos::as; using Teuchos::incrVerbLevel; typedef ScalarTraits<Scalar> ST; typedef Thyra::NonlinearSolverBase<Scalar> NSB; typedef Teuchos::VerboseObjectTempState<NSB> VOTSNSB; RCP<FancyOStream> out = this->getOStream(); Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel(); Teuchos::OSTab ostab(out,1,"takeStep"); VOTSNSB solver_outputTempState(solver_,out,incrVerbLevel(verbLevel,-1)); if ( !is_null(out) && as<int>(verbLevel) >= as<int>(Teuchos::VERB_LOW) ) { *out << "\nEntering " << Teuchos::TypeNameTraits<ImplicitRKStepper<Scalar> >::name() << "::takeStep("<<dt<<","<<toString(stepSizeType)<<") ...\n"; } if (!isInitialized_) { initialize_(); } TEUCHOS_TEST_FOR_EXCEPT( stepSizeType != STEP_TYPE_FIXED ); // ToDo: Handle variable case later // A) Set up the IRK ModelEvaluator so that it can represent the time step // equation to be solved. // Set irkModel_ with x_old_, t_old_, and dt V_V( x_old_.ptr(), *x_ ); Scalar current_dt = dt; Scalar t = timeRange_.upper(); // B) Solve the timestep equation // Set the guess for the stage derivatives to zero (unless we can think of // something better) V_S( Teuchos::rcp_dynamic_cast<Thyra::VectorBase<Scalar> >(x_stage_bar_).ptr(), ST::zero() ); if (!isDirk_) { // General Implicit RK Case: RCP<ImplicitRKModelEvaluator<Scalar> > firkModel_ = Teuchos::rcp_dynamic_cast<ImplicitRKModelEvaluator<Scalar> >(irkModel_,true); firkModel_->setTimeStepPoint( x_old_, t, current_dt ); // Solve timestep equation solver_->solve( &*x_stage_bar_ ); } else { // Diagonal Implicit RK Case: RCP<DiagonalImplicitRKModelEvaluator<Scalar> > dirkModel_ = Teuchos::rcp_dynamic_cast<DiagonalImplicitRKModelEvaluator<Scalar> >(irkModel_,true); dirkModel_->setTimeStepPoint( x_old_, t, current_dt ); int numStages = irkButcherTableau_->numStages(); for (int stage=0 ; stage < numStages ; ++stage) { dirkModel_->setCurrentStage(stage); solver_->solve( &*(x_stage_bar_->getNonconstVectorBlock(stage)) ); dirkModel_->setStageSolution( stage, *(x_stage_bar_->getVectorBlock(stage)) ); } } // C) Complete the step ... // Combine the stage derivatives with the Butcher tableau "b" vector to obtain the solution at the final time. // x_{k+1} = x_k + dt*sum_{i}^{p}(b_i*x_stage_bar_[i]) assembleIRKSolution( irkButcherTableau_->b(), current_dt, *x_old_, *x_stage_bar_, outArg(*x_) ); // Update time range timeRange_ = timeRange(t,t+current_dt); numSteps_++; return current_dt; }
void DefaultIntegrator<Scalar>::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 ) { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR_DIFF("Rythmos:DefaultIntegrator::getFwdPoints", TopLevel); #endif using Teuchos::incrVerbLevel; #ifndef _MSC_VER using Teuchos::Describable; #endif typedef Teuchos::ScalarTraits<Scalar> ST; typedef InterpolationBufferBase<Scalar> IBB; typedef Teuchos::VerboseObjectTempState<IBB> VOTSIBB; finalizeSetup(); RCP<Teuchos::FancyOStream> out = this->getOStream(); Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel(); Teuchos::OSTab tab(out); VOTSIBB stepper_outputTempState(stepper_,out,incrVerbLevel(verbLevel,-1)); if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nEntering " << this->Describable::description() << "::getFwdPoints(...) ...\n" << "\nStepper: " << Teuchos::describe(*stepper_,verbLevel); if ( includesVerbLevel(verbLevel,Teuchos::VERB_MEDIUM) ) *out << "\nRequested time points: " << Teuchos::toString(time_vec) << "\n"; // // 0) Initial setup // const int numTimePoints = time_vec.size(); // Assert preconditions assertTimePointsAreSorted(time_vec); TEST_FOR_EXCEPT(accuracy_vec!=0); // ToDo: Remove accuracy_vec! // Resize the storage for the output arrays if (x_vec) x_vec->resize(numTimePoints); if (xdot_vec) xdot_vec->resize(numTimePoints); // This int records the next time point offset in time_vec[timePointIndex] // that needs to be handled. This gets updated as the time points are // filled below. int nextTimePointIndex = 0; assertNoTimePointsBeforeCurrentTimeRange(*this,time_vec,nextTimePointIndex); // // 1) First, get all time points that fall within the current time range // { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: getPoints"); #endif // 2007/10/05: rabartl: ToDo: Get points from trailingInterpBuffer_ first! getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex); } // // 2) Advance the stepper to satisfy time points in time_vec that fall // before the current time. // while ( nextTimePointIndex < numTimePoints ) { // Use the time stepping algorithm to step up to or past the next // requested time but not so far as to step past the point entirely. const Scalar t = time_vec[nextTimePointIndex]; bool advanceStepperToTimeSucceeded = false; { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: advanceStepperToTime"); #endif advanceStepperToTimeSucceeded= advanceStepperToTime(t); } if (!advanceStepperToTimeSucceeded) { bool reachedMaxNumTimeSteps = (currTimeStepIndex_ >= maxNumTimeSteps_); if (reachedMaxNumTimeSteps) { // Break out of the while loop and attempt to exit gracefully. break; } TEST_FOR_EXCEPTION( !advanceStepperToTimeSucceeded, Exceptions::GetFwdPointsFailed, this->description() << "\n\n" "Error: The integration failed to get to time " << t << " and only achieved\n" "getting to " << stepper_->getTimeRange().upper() << "!" ); } // Extract the next set of points (perhaps just one) from the stepper { #ifdef ENABLE_RYTHMOS_TIMERS TEUCHOS_FUNC_TIME_MONITOR("Rythmos:DefaultIntegrator::getFwdPoints: getPoints (fwd)"); #endif getCurrentPoints(*stepper_,time_vec,x_vec,xdot_vec,&nextTimePointIndex); } } if ( includesVerbLevel(verbLevel,Teuchos::VERB_LOW) ) *out << "\nLeaving " << this->Describable::description() << "::getFwdPoints(...) ...\n"; }