2#ifndef RYTHMOS_INTEGRATION_OBSERVER_BASE_HPP
3#define RYTHMOS_INTEGRATION_OBSERVER_BASE_HPP
5#include "Rythmos_Types.hpp"
6#include "Teuchos_Describable.hpp"
7#include "Teuchos_VerboseObject.hpp"
13template<
class Scalar>
class TimeRange;
14template<
class Scalar>
class StepControlInfo;
15template<
class Scalar>
class StepperBase;
25 :
virtual public Teuchos::Describable,
26 virtual public Teuchos::VerboseObject<IntegrationObserverBase<Scalar> >
35 virtual RCP<IntegrationObserverBase<Scalar> >
116 const int timeStepIter
153 const int timeStepIter
194 const int timeStepIter
205template<
class Scalar>
206bool isInitialTimeStep(
211 typedef Teuchos::ScalarTraits<Scalar> ST;
212 return compareTimeValues(currentTimeRange.
lower(), fullTimeRange.
lower()) == ST::zero();
221template<
class Scalar>
223 const TimeRange<Scalar> ¤tTimeRange,
224 const TimeRange<Scalar> &fullTimeRange
227 typedef Teuchos::ScalarTraits<Scalar> ST;
228 return compareTimeValues(currentTimeRange.upper(), fullTimeRange.upper()) >= ST::zero();
242template<
class Scalar>
249template<
class Scalar>
256template<
class Scalar>
267template<
class Scalar>
Base class for strategy objects that observe and time integration by observing the stepper object.
virtual void observeCompletedTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)=0
Observe a successfully completed integration step.
virtual void observeStartTimeIntegration(const StepperBase< Scalar > &stepper)
Observe the beginning of a time integration loop.
virtual void observeEndTimeIntegration(const StepperBase< Scalar > &stepper)
Observe the end of a time integration loop.
virtual void resetIntegrationObserver(const TimeRange< Scalar > &integrationTimeDomain)=0
Reset the observer to prepair it to observe another integration.
virtual void observeStartTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)
Observer the beginning of an integration step.
virtual void observeFailedTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)
Observer a failed integration step.
virtual RCP< IntegrationObserverBase< Scalar > > cloneIntegrationObserver() const =0
Clone this integration observer if supported .
Simple struct to aggregate integration/stepper control information.
Base class for defining stepper functionality.