29 #ifndef Rythmos_LOGGING_INTEGRATION_OBSERVER_HPP 
   30 #define Rythmos_LOGGING_INTEGRATION_OBSERVER_HPP 
   32 #include "Rythmos_IntegrationObserverBase.hpp" 
   33 #include "Teuchos_RCP.hpp" 
   43 template<
class Scalar>
 
   50   void resetLogCounters();
 
   52   Teuchos::RCP<const std::map<std::string,int> > getCounters();
 
   54   Teuchos::RCP<const std::list<std::string> > getOrder();
 
   71     const int timeStepIter
 
   77     const int timeStepIter
 
   83     const int timeStepIter
 
   94   const std::string nameCloneIntegrationObserver_;
 
   95   const std::string nameResetIntegrationObserver_;
 
   96   const std::string nameObserveStartTimeIntegration_;
 
   97   const std::string nameObserveEndTimeIntegration_;
 
   98   const std::string nameObserveStartTimeStep_;
 
   99   const std::string nameObserveCompletedTimeStep_;
 
  100   const std::string nameObserveFailedTimeStep_;
 
  112   void logCall(
const std::string call) 
const;
 
  116   Teuchos::RCP< std::map<std::string,int> > counters_;
 
  117   Teuchos::RCP< std::list<std::string> > order_;
 
  126 template<
class Scalar>
 
  127 Teuchos::RCP<LoggingIntegrationObserver<Scalar> >
 
  130   const Teuchos::RCP<LoggingIntegrationObserver<Scalar> > lio = 
 
  140 template<
typename Scalar>
 
  142   nameCloneIntegrationObserver_(
"cloneIntegrationObserver"),
 
  143   nameResetIntegrationObserver_(
"resetIntegrationObserver"),
 
  144   nameObserveStartTimeIntegration_(
"observeStartTimeIntegration"),
 
  145   nameObserveEndTimeIntegration_(
"observeEndTimeIntegration"),
 
  146   nameObserveStartTimeStep_(
"observeStartTimeStep"),
 
  147   nameObserveCompletedTimeStep_(
"observeCompletedTimeStep"),
 
  148   nameObserveFailedTimeStep_(
"observeFailedTimeStep")
 
  150   counters_ = Teuchos::rcp(
new std::map<std::string,int>);
 
  151   order_ = Teuchos::rcp(
new std::list<std::string>);
 
  152   this->resetLogCounters();
 
  155 template<
typename Scalar>
 
  156 void LoggingIntegrationObserver<Scalar>::
 
  159   (*counters_)[nameCloneIntegrationObserver_] = 0;
 
  160   (*counters_)[nameResetIntegrationObserver_] = 0;
 
  161   (*counters_)[nameObserveStartTimeIntegration_] = 0;
 
  162   (*counters_)[nameObserveEndTimeIntegration_] = 0;
 
  163   (*counters_)[nameObserveStartTimeStep_] = 0;
 
  164   (*counters_)[nameObserveCompletedTimeStep_] = 0;
 
  165   (*counters_)[nameObserveFailedTimeStep_] = 0;
 
  169 template<
typename Scalar>
 
  170 RCP<IntegrationObserverBase<Scalar> > 
 
  173   logCall(nameCloneIntegrationObserver_);
 
  174   Teuchos::RCP<IntegrationObserverBase<Scalar> > observer = 
 
  179 template<
typename Scalar>
 
  184   logCall(nameResetIntegrationObserver_);
 
  187 template<
typename Scalar>
 
  191   logCall(nameObserveStartTimeIntegration_);
 
  194 template<
typename Scalar>
 
  198   logCall(nameObserveEndTimeIntegration_);
 
  201 template<
typename Scalar>
 
  205     const int timeStepIter
 
  208   logCall(nameObserveStartTimeStep_);
 
  211 template<
typename Scalar>
 
  215     const int timeStepIter
 
  218   logCall(nameObserveCompletedTimeStep_);
 
  221 template<
typename Scalar>
 
  225     const int timeStepIter
 
  228   logCall(nameObserveFailedTimeStep_);
 
  231 template<
typename Scalar>
 
  237 template<
typename Scalar>
 
  238 RCP<const std::list<std::string> > LoggingIntegrationObserver<Scalar>::getOrder()
 
  244 template<
typename Scalar>
 
  245 void LoggingIntegrationObserver<Scalar>::logCall(
const std::string call)
 const 
  247   (*counters_)[call] += 1;
 
  248   order_->push_back(call);
 
  255 #endif //Rythmos_LOGGING_INTEGRATION_OBSERVER_HPP 
void observeStartTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)
Observer the beginning of an integration step. 
 
Base class for defining stepper functionality. 
 
Simple struct to aggregate integration/stepper control information. 
 
void observeEndTimeIntegration(const StepperBase< Scalar > &stepper)
Observe the end of a time integration loop. 
 
void resetIntegrationObserver(const TimeRange< Scalar > &integrationTimeDomain)
Reset the observer to prepair it to observe another integration. 
 
RCP< IntegrationObserverBase< Scalar > > cloneIntegrationObserver() const 
Clone this integration observer if supported . 
 
Base class for strategy objects that observe and time integration by observing the stepper object...
 
Logging IntegrationOberserver that counts calls to observer functions and lists their order...
 
void observeFailedTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)
Observer a failed integration step. 
 
Teuchos::RCP< LoggingIntegrationObserver< Scalar > > createLoggingIntegrationObserver()
Nonmember constructor. 
 
void observeStartTimeIntegration(const StepperBase< Scalar > &stepper)
Observe the beginning of a time integration loop. 
 
void observeCompletedTimeStep(const StepperBase< Scalar > &stepper, const StepControlInfo< Scalar > &stepCtrlInfo, const int timeStepIter)
Observe a successfully completed integration step.