MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
IterationPack_Algorithm.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ALGORITHM_H
43 #define ALGORITHM_H
44 
45 #include <assert.h>
46 
47 #include <string>
48 #include <deque>
49 #include <list>
50 #include <vector>
51 #include <algorithm>
52 
56 #include "Teuchos_RCP.hpp"
57 #include "Teuchos_Assert.hpp"
59 
60 namespace IterationPack {
61 
62 // ToDo: 7/31/98: Finish documentation.
63 
114 class Algorithm {
115 public:
116 
119 
127  typedef size_t poss_type;
129  enum { DOES_NOT_EXIST = 1000 }; // never be that many steps
132 
134  class DoesNotExist : public std::logic_error
135  {public: DoesNotExist(const std::string& what_arg) : std::logic_error(what_arg) {}};
136 
138  class AlreadyExists : public std::logic_error
139  {public: AlreadyExists(const std::string& what_arg) : std::logic_error(what_arg) {}};
140 
142  class InvalidControlProtocal : public std::logic_error
143  {public: InvalidControlProtocal(const std::string& what_arg) : std::logic_error(what_arg) {}};
144 
146  class InvalidRunningState : public std::logic_error
147  {public: InvalidRunningState(const std::string& what_arg) : std::logic_error(what_arg) {}};
148 
150  class InvalidConfigChange : public std::logic_error
151  {public: InvalidConfigChange(const std::string& what_arg) : std::logic_error(what_arg) {}};
152 
154  class AlgorithmInterrupted : public std::runtime_error
155  {public: AlgorithmInterrupted(const std::string& what_arg) : std::runtime_error(what_arg) {}};
156 
158 
161 
165  Algorithm();
166 
168  virtual ~Algorithm();
169 
171 
206  STANDARD_MEMBER_COMPOSITION_MEMBERS( std::string, interrupt_file_name );
207 
210 
212  void set_state(const state_ptr_t& state);
216  const state_ptr_t& get_state() const;
220  const AlgorithmState& state() const;
221 
223 
226 
228  void set_track(const track_ptr_t& track);
232  const track_ptr_t& get_track() const;
236  const AlgorithmTracker& track() const;
237 
239 
242 
244  virtual void max_iter(size_t max_iter);
246  virtual size_t max_iter() const;
247 
249 
252 
257  virtual void max_run_time(double max_iter);
259  virtual double max_run_time() const;
260 
262 
265 
267  virtual int num_steps() const;
268 
274  virtual poss_type get_step_poss(const std::string& step_name) const;
275 
282  virtual const std::string& get_step_name(poss_type step_poss) const;
283 
290  virtual step_ptr_t& get_step(poss_type step_poss);
291 
293  virtual const step_ptr_t& get_step(poss_type step_poss) const;
294 
296 
299 
306  virtual int num_assoc_steps(poss_type step_poss, EAssocStepType type) const;
307 
317  virtual poss_type get_assoc_step_poss(poss_type step_poss, EAssocStepType type
318  ,const std::string& assoc_step_name) const;
319 
328  virtual const std::string& get_assoc_step_name(poss_type step_poss, EAssocStepType type
329  , poss_type assoc_step_poss) const;
330 
340  virtual step_ptr_t& get_assoc_step(poss_type step_poss, EAssocStepType type
341  , poss_type assoc_step_poss);
342 
344  virtual const step_ptr_t& get_assoc_step(poss_type step_poss, EAssocStepType type
345  , poss_type assoc_step_poss) const;
346 
348 
351 
364  virtual void insert_step(poss_type step_poss, const std::string& step_name, const step_ptr_t& step);
365 
375  virtual void change_step_name(poss_type step_poss, const std::string& new_name);
376 
386  virtual void replace_step(poss_type step_poss, const step_ptr_t& step);
387 
398  virtual void remove_step(poss_type step_poss);
399 
401 
404 
419  virtual void insert_assoc_step(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss
420  , const std::string& assoc_step_name, const step_ptr_t& assoc_step);
421 
435  virtual void remove_assoc_step(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss);
436 
438 
441 
444 
457  virtual void begin_config_update();
458 
471  virtual void end_config_update();
472 
474 
477 
486  virtual void do_step_next(const std::string& step_name);
487 
496  virtual void do_step_next(poss_type step_poss);
497 
504  virtual const std::string& what_is_next_step_name() const;
505 
512  virtual poss_type what_is_next_step_poss() const;
513 
530  virtual bool do_step(const std::string& step_name);
531 
547  virtual bool do_step(poss_type step_poss);
548 
558  virtual void terminate(bool success);
559 
561 
564 
598  virtual EAlgoReturn do_algorithm(poss_type step_poss = 1);
599 
601 
604 
607  virtual void print_steps(std::ostream& out) const;
608 
611  virtual void print_algorithm(std::ostream& out) const;
612 
614 
618 
625  virtual void set_algo_timing( bool algo_timing );
626 
628  virtual bool algo_timing() const;
629 
638  virtual void print_algorithm_times( std::ostream& out ) const;
639 
655  void get_step_times_k( int offset, double step_times[] ) const;
656 
660  void get_final_step_stats( size_t step, double* total, double* average, double* min, double* max, double* percent) const;
661 
663 
664 private:
665 
666  // /////////////////////////////////////////////////////
667  // Private types
668 
670  template<class T_Step_ptr>
673  step_ptr_and_name(const T_Step_ptr& _step_ptr
674  , const std::string& _name )
675  : step_ptr(_step_ptr), name(_name)
676  {}
678  T_Step_ptr step_ptr;
679  //
680  std::string name;
681  }; // end struct step_ptr_and_name
682 
686  typedef std::deque<steps_ele_t> steps_t;
687 
691  typedef std::list<assoc_steps_ele_list_ele_t> assoc_steps_ele_list_t;
696  { return assoc_steps_lists_[i]; }
698  const assoc_steps_ele_list_t& operator[](int i) const
699  { return assoc_steps_lists_[i]; }
700  private:
702  };
703 
704  //typedef assoc_steps_ele_list_t[2] assoc_steps_ele_t; // PRE_STEP, POST_STEP
706  typedef std::deque<assoc_steps_ele_t> assoc_steps_t;
707 
710 
712  template<class T_ele>
713  class name_comp {
714  public:
716  name_comp(const std::string& name) : name_(name) {}
718  bool operator()(const T_ele& ele) { return ele.name == name_; }
719  private:
720  const std::string& name_;
721  }; // end class name_comp
722 
723  typedef std::vector<double> step_times_t;
724 
726  static const int
733  enum { NUM_STEP_TIME_STATS = 5 };
734 
737  Algorithm& algo
738  ,poss_type step_poss
739  ,EDoStepType type
740  ,poss_type assoc_step_poss
741  );
742 
743  // /////////////////////////////////////////////////////
744  // Private data members
745 
746  // aggregate members
747 
748 #ifdef DOXYGEN_COMPILE
752 #else
753  state_ptr_t state_;
754  // RCP<...> object for the aggragate AlgorithmState object.
755 
756  track_ptr_t track_;
757  // RCP<...> object for the aggragate AlgorithmTracker object.
758 #endif
759 
760  // algorithm control etc.
761 
763  // The state this Algorithm object is in:
764  //
765  // NOT_RUNNING do_algorithm() is not active.
766  // RUNNING do_algorithm() has been called and is active.
767  // RUNNING_BEING_CONFIGURED
768  // do_algorithm() is active and begin_config_update() has been called
769  // but end_config_update() has not.
770  //
771  // Note: Only change this variable through the private function change_running_state(...)
772  // unless you are 100% sure that you know what you are doing!
773 
774  size_t first_k_;
775  // The first iteration from state().k().
776 
777  size_t max_iter_;
778  // The maximum number of iterations that <tt>this</tt> will execute.
779 
781  // The maximum amount of time the algorithm is allowed to execute.
782 
784  // Flag for if it is time to terminate do_algorithm().
785 
787  // The next step possition that <tt>this</tt> will execute when control is returned to do_algorithm().
788 
789  const std::string* next_step_name_;
790  // The name of the next step that <tt>this</tt> will execute when control is returned to do_algorithm().
791 
793  // Flag for if do_step_next() was called so that <tt>do_algorithm()</tt> can validate
794  // that if a step object returned <tt>false</tt> from its <tt>do_step()</tt> operation that it
795  // must also call <tt>do_step_next()</tt> to specify a step to jump to.
796 
798  // The current step being executed in do_algorithm().
799  // If the current step being executed is changed during the imp_do_step() operation, then
800  // imp_do_step() must adjust to this step.
801 
803  // The name of the current step that is saved when begin_config_update() is called
804  // so that curr_step_poss_ can be reset when end_config_update() is called.
805 
807  // The name of the next step to call so that when begin_config_update() is called
808  // so that next_step_poss_ and next_step_name_ can be reset when end_config_update()
809  // is called.
810 
812  // A flag that is set to true when a runtime configuration has been preformed. It
813  // is used to flag this event for imp_do_assoc_steps().
814 
815  // step and associated step object data structures
816 
818  // Array of std::pair<RCP<step_ptr_t>,std::string> objects.
819  //
820  // *steps_[step_poss].first returns the step object for step_poss = 1...steps_.size().
821  // steps_[step_poss].second returns the name of the step for step_poss = 1...steps_.size().
822 
824  // Array of two lists of std::pair<step_ptr_t,std::string> objects
825  //
826  // *(assoc_steps_[step_poss][PRE_STEP].begin() + pre_step_poss).first gives the pre step object.
827  // (assoc_steps_[step_poss][PRE_STEP].begin() + pre_step_poss).second gives the name of the pre step
828  // *(assoc_steps_[step_poss][POST_STEP].begin() + post_step_poss).first gives the post step object.
829  // (assoc_steps_[step_poss][POST_STEP].begin() + post_step_poss).second gives the name of the post step
830 
832  // If true each step will be timed.
833 
835  // Array of step times ( size (max_iter() + 1 + NUM_STEP_TIME_STATS) * (num_steps() + 1) ).
836  // The time in sec. for step step_i (one based)
837  // for iteration iter_k (zero based) is:
838  // step_times_[ iter_k + (step_i - 1) * (max_iter() + 1 + NUM_STEP_TIME_STATS) ].
839  // Therefore the times for each step are stored by column (consecutive elements)
840  // so that statistics will be easy to compute at the end.
841  // The last five elements after max_iter() for each step are reserved for:
842  // * total time for the step
843  // * average time for the step
844  // * min step time
845  // * max step time
846  // * percentage for each step to the total.
847  // The last column is for the total times for each iteration with the last five
848  // elements being for the statistics for each iteration.
849 
850  mutable bool time_stats_computed_;
851  // A flag for if the timing statistics have already been computed or not.
852 
853  mutable double total_time_;
854  // Records the total computed time for the algorithm.
855 
856  // /////////////////////////////////////////////////////
857  // Private member functions
858 
860  poss_type validate(poss_type step_poss, int past_end = 0) const;
861 
863  poss_type validate(const assoc_steps_ele_list_t& assoc_list, poss_type assoc_step_poss, int past_end = 0) const;
864 
867 
870 
873 
875  void validate_not_curr_step(poss_type step_poss) const;
876 
878  void validate_not_next_step(const std::string& step_name) const;
879 
882  steps_t::iterator step_itr(const std::string& step_name);
883 
885  steps_t::const_iterator step_itr(const std::string& step_name) const;
886 
889  steps_t::iterator step_itr_and_assert(const std::string& step_name);
890 
892  steps_t::const_iterator step_itr_and_assert(const std::string& step_name) const;
893 
896  static assoc_steps_ele_list_t::iterator assoc_step_itr(assoc_steps_ele_list_t& assoc_list
897  , const std::string& assoc_step_name);
898 
900  static assoc_steps_ele_list_t::const_iterator assoc_step_itr(const assoc_steps_ele_list_t& assoc_list
901  , const std::string& assoc_step_name);
902 
904  bool imp_do_step(poss_type step_poss);
905 
908 
910  void imp_inform_steps(inform_func_ptr_t inform_func_ptr);
911 
913  void imp_print_algorithm(std::ostream& out, bool print_steps) const;
914 
916  EDoStepType do_step_type(EAssocStepType assoc_step_type);
917 
920 
922  void compute_final_time_stats() const;
923 
924  // Look for interrup
925  void look_for_interrupt();
926 
927 public:
928 
929  // This is put here out of need. Not for any user to call!
930  static void interrupt();
931 
932 }; // end class Algorithm
933 
934 // //////////////////////////////////////////////////////////////////////////////////////////////////
935 // Inline member function definitions for Algorithm
936 
937 // «std comp» members for state
938 
939 inline
941 { state_ = state; }
942 
943 inline
945 { return state_; }
946 
947 inline
949 { return state_; }
950 
951 inline
953 { TEUCHOS_TEST_FOR_EXCEPT( !( state_.get() ) ); return *state_; }
954 
955 inline
956 const AlgorithmState& Algorithm::state() const
957 { TEUCHOS_TEST_FOR_EXCEPT( !( state_.get() ) ); return *state_; }
958 
959 // «std comp» members for track
960 
961 inline
963 { track_ = track; }
964 
965 inline
967 { return track_; }
968 
969 inline
971 { return track_; }
972 
973 inline
975 { TEUCHOS_TEST_FOR_EXCEPT( !( track_.get() ) ); return *track_; }
976 
977 inline
978 const AlgorithmTracker& Algorithm::track() const
979 { TEUCHOS_TEST_FOR_EXCEPT( !( track_.get() ) ); return *track_; }
980 
981 // running state
982 
983 inline
985 { return running_state_; }
986 
987 // lookup iterator given name
988 
989 inline
990 Algorithm::steps_t::iterator Algorithm::step_itr(const std::string& step_name)
991 {
992  return std::find_if( steps_.begin() , steps_.end()
993  , name_comp<steps_ele_t>(step_name) );
994 }
995 
996 inline
997 Algorithm::steps_t::const_iterator Algorithm::step_itr(const std::string& step_name) const
998 {
999  return std::find_if( steps_.begin() , steps_.end()
1000  , name_comp<steps_ele_t>(step_name) );
1001 }
1002 
1003 inline
1004 Algorithm::assoc_steps_ele_list_t::iterator Algorithm::assoc_step_itr(
1005  assoc_steps_ele_list_t& assoc_list, const std::string& assoc_step_name)
1006 {
1007  return std::find_if( assoc_list.begin() , assoc_list.end()
1008  , name_comp<assoc_steps_ele_list_ele_t>(assoc_step_name) );
1009 }
1010 
1011 inline
1012 Algorithm::assoc_steps_ele_list_t::const_iterator Algorithm::assoc_step_itr(
1013  const assoc_steps_ele_list_t& assoc_list, const std::string& assoc_step_name)
1014 {
1015  return std::find_if( assoc_list.begin() , assoc_list.end()
1016  , name_comp<assoc_steps_ele_list_ele_t>(assoc_step_name) );
1017 }
1018 
1019 inline
1021  switch(assoc_step_type) {
1022  case PRE_STEP : return DO_PRE_STEP;
1023  case POST_STEP : return DO_POST_STEP;
1024  }
1025  TEUCHOS_TEST_FOR_EXCEPT( !( true ) );
1026  return DO_PRE_STEP; // will never execute.
1027 }
1028 
1029 } // end namespace IterationPack
1030 
1031 #endif // ALGORITHM_H
Base type for all objects that perform steps in an Algorithm.
std::vector< double > step_times_t
void validate_not_in_state(ERunningState running_state) const
Validate that this is not in a specific running state.
virtual EAlgoReturn do_algorithm(poss_type step_poss=1)
Called by clients to begin an algorithm.
step_ptr_and_name< step_ptr_t > assoc_steps_ele_list_ele_t
step_ptr_and_name< step_ptr_t > steps_ele_t
Teuchos::RCP< AlgorithmTracker > track_ptr_t
steps_t::iterator step_itr(const std::string &step_name)
Find a step given its name and throw a DoesNotExist exception if not found.
bool imp_do_step(poss_type step_poss)
Algorithm()
Constructs an algorithm with no steps and a default of max_iter() == 100.
virtual poss_type get_step_poss(const std::string &step_name) const
Return the possition in the major loop of a named step.
ERunningState running_state() const
Return the current running state of this algorithm object.
virtual void print_steps(std::ostream &out) const
Print out just a listing of the steps, their positions in the algorithm and the subclasses.
virtual void begin_config_update()
Changes from running_state() == RUNNING to running_state() == RUNNING_BEING_CONFIGURED.
virtual void insert_assoc_step(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss, const std::string &assoc_step_name, const step_ptr_t &assoc_step)
Insert an pre or post step into for the main step step_poss into the possition assoc_step_poss.
Thrown if name or id does not exist.
typedef void(CALL_API *RTOp_reduct_op_func_ptr_t)(void *
Thrown if an invalid control protocal is used.
void get_step_times_k(int offset, double step_times[]) const
Returns the step_times for iteration offset.
virtual const std::string & what_is_next_step_name() const
Returns the name of the next step this will call the next time it calls a step.
void get_final_step_stats(size_t step, double *total, double *average, double *min, double *max, double *percent) const
Returns the final statistics for a given step Do not call when algorithm is running.
virtual void insert_step(poss_type step_poss, const std::string &step_name, const step_ptr_t &step)
Insert a step object with the name step_name into the possition step_poss.
void imp_inform_steps(inform_func_ptr_t inform_func_ptr)
virtual int num_steps() const
Return the number of main steps.
void validate_in_state(ERunningState running_state) const
Validate that this is in a specific running state.
virtual void change_step_name(poss_type step_poss, const std::string &new_name)
Change the name of an existing step.
STANDARD_MEMBER_COMPOSITION_MEMBERS(std::string, interrupt_file_name)
Name of an file that will cause the algorithm to terminate.
std::deque< assoc_steps_ele_t > assoc_steps_t
void validate_not_next_step(const std::string &step_name) const
Validate that the step_name in not the next step.
virtual void print_algorithm(std::ostream &out) const
Print out the entire algorithm by calling print_step(...) on the step objects.
virtual bool do_step(const std::string &step_name)
Calls do_step() on all of the pre step objects the step object and the post step objects in order for...
virtual void terminate(bool success)
Called by step objects to terminate the algorithm.
void set_state(const state_ptr_t &state)
void imp_print_algorithm(std::ostream &out, bool print_steps) const
const assoc_steps_ele_list_t & operator[](int i) const
virtual void replace_step(poss_type step_poss, const step_ptr_t &step)
Replace the step object of an existing step.
std::ostream * out
steps_t::iterator step_itr_and_assert(const std::string &step_name)
Find a step given its name and throw a DoesNotExist exception if not found.
virtual void do_step_next(const std::string &step_name)
Called by step objects to set the step (given its name) that this will envoke the next time this call...
virtual poss_type what_is_next_step_poss() const
Returns the possition of the next step this will call the next time it calls a step.
Abstacts a set of iteration quantities for an iterative algorithm.
virtual poss_type get_assoc_step_poss(poss_type step_poss, EAssocStepType type, const std::string &assoc_step_name) const
Return the possition of the pre or post step for the main step_poss.
AlgorithmTracker & track()
virtual int num_assoc_steps(poss_type step_poss, EAssocStepType type) const
Return the number of pre or post steps for the main step step_poss.
Teuchos::RCP< AlgorithmState > state_ptr_t
virtual step_ptr_t & get_assoc_step(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss)
Return the RCP<...> object for the associated step object at step_poss and assoc_step_poss.
virtual step_ptr_t & get_step(poss_type step_poss)
Return the RCP<...> object for the step object at step_poss.
std::list< assoc_steps_ele_list_ele_t > assoc_steps_ele_list_t
AlgorithmState & state()
static assoc_steps_ele_list_t::iterator assoc_step_itr(assoc_steps_ele_list_t &assoc_list, const std::string &assoc_step_name)
Find a an associated step given its name and throw a DoesNotExist exception if not found...
void validate_not_curr_step(poss_type step_poss) const
Validate that the step_poss in not the current step.
EDoStepType do_step_type(EAssocStepType assoc_step_type)
EAssocStepType -> EDoStepType.
virtual void end_config_update()
Changes from running_state() == RUNNING_BEING_CONFIGURED to running_state() == RUNNING.
std::deque< steps_ele_t > steps_t
void(AlgorithmStep::* inform_func_ptr_t)(Algorithm &algo, poss_type step_poss, EDoStepType type, poss_type assoc_step_poss)
virtual const std::string & get_step_name(poss_type step_poss) const
Return the name of a step given its possition.
virtual double max_run_time() const
EAlgoReturn finalize_algorithm(EAlgoReturn algo_return)
Thrown if a member function is called while this is in an invalid running state.. ...
Acts as the central hub for an iterative algorithm.
Thrown if Algorithm was interrupted by the user.
step_ptr_and_name(const T_Step_ptr &_step_ptr, const std::string &_name)
virtual void remove_assoc_step(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss)
Remove an pre or post step for the main step step_poss in the possition assoc_step_poss.
virtual void print_algorithm_times(std::ostream &out) const
Outputs table of times for each step, cummulative times and other statistics.
Used to ouput iteration results and other information.
virtual const std::string & get_assoc_step_name(poss_type step_poss, EAssocStepType type, poss_type assoc_step_poss) const
Return the name of the pre or post step at step_poss and at assoc_step_poss.
Thrown if a member function is called while this is in an invalid running state.. ...
virtual void remove_step(poss_type step_poss)
Remove an existing step object and all of its pre and post steps.
bool imp_do_assoc_steps(EAssocStepType type)
virtual size_t max_iter() const
Teuchos::RCP< AlgorithmStep > step_ptr_t
virtual void set_algo_timing(bool algo_timing)
Causes algorithm to be timed.
void set_track(const track_ptr_t &track)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
void change_running_state(ERunningState running_state)
Change the running state.
poss_type validate(poss_type step_poss, int past_end=0) const
Validate a step_poss and throw a DoesNotExist exception if it does not.