DOpE
spacetimehandler_base.h
Go to the documentation of this file.
1 
24 #ifndef _SPACE_TIME_HANDLER_BASE_H_
25 #define _SPACE_TIME_HANDLER_BASE_H_
26 
27 #include <lac/vector.h>
28 #include <lac/block_vector_base.h>
29 #include <lac/block_vector.h>
30 
31 #include <vector>
32 #include <iostream>
33 #include <sstream>
34 #include <limits>
35 
36 #include "timedofhandler.h"
37 #include "timeiterator.h"
38 #include "dopetypes.h"
39 
40 namespace DOpE
41 {
46 template<typename VECTOR>//has to be a class template because you can not have virtual member function templates.
48 {
49  public:
50 
51  SpaceTimeHandlerBase(DOpEtypes::ControlType control_type = DOpEtypes::undefined) : _control_type(control_type)
52  {
53  _state_ticket = 1;
54  _control_ticket = 1;
55  }
56 
57  SpaceTimeHandlerBase(const dealii::Triangulation<1> & times, DOpEtypes::ControlType type = DOpEtypes::undefined) :
58  _tdfh(times), _interval(_tdfh.first_interval()), _control_type(type)
59  {
60  _state_ticket = 1;
61  _control_ticket = 1;
62  }
63 
64  SpaceTimeHandlerBase(const dealii::Triangulation<1> & times,
65  const dealii::FiniteElement<1>& fe,
67  _tdfh(times, fe), _interval(_tdfh.first_interval()), _control_type(type)
68  {
69  _state_ticket = 1;
70  _control_ticket = 1;
71  }
72 
73 
75  {
76  _tdfh.clear();
77  }
78 
82  void ReInitTime()
83  {
85  _tdfh.distribute_dofs();
86  //FIXME When we have temporal discretization also for control and constraint,
87  //one has to increment here the _control_ticket!
88  }
89 
97  unsigned int GetMaxTimePoint() const
98  {
99  return (_tdfh.GetNbrOfDoFs()-1);//because we start counting at 0.
100  }
101 
105  unsigned int GetNbrOfIntervals() const
106  {
107  return _tdfh.GetNbrOfIntervals();
108  }
109 
110 
116  void SetInterval(const TimeIterator& it)
117  {
118  _interval = it;
119  }
120 
121 
127  const TimeIterator& GetInterval() const
128  {
129  return _interval;
130  }
131 
132 
139  double GetTime(unsigned int time_point) const
140  {
141 
142  return _tdfh.GetTime(time_point);
143  }
144 
145 
150  {
151  return _tdfh;
152  }
153 
154 
158  const std::vector<double>& GetTimes() const
159  {
160  return _tdfh.GetTimes();
161  }
162 
174  void
175  GetTimes(const TimeIterator& interval, std::vector<double>& local_times) const
176  {
177  return _tdfh.GetTimes(interval, local_times);
178  }
179 
180 
191  bool IsValidStateTicket(unsigned int& ticket) const
192  {
193  bool ret = (ticket == _state_ticket);
194  ticket = _state_ticket;
195  return ret;
196  }
197 
208  bool IsValidControlTicket(unsigned int& ticket) const
209  {
210  bool ret = (ticket == _control_ticket);
211  ticket = _control_ticket;
212  return ret;
213  }
214 
215 
221  {
222  return _control_type;
223  }
224 
225 
237  virtual void
238  InterpolateControl(VECTOR& /*result*/, const std::vector<VECTOR*> &/*local_vectors*/,
239  double /*t*/, const TimeIterator&/*interval*/) const { abort(); }
251  virtual void
252  InterpolateState(VECTOR& /*result*/, const std::vector<VECTOR*> &/*local_vectors*/,
253  double /*t*/, const TimeIterator&/*interval*/) const = 0;
265  virtual void InterpolateConstraint(VECTOR& /*result*/, const std::vector<VECTOR*> &/*local_vectors*/,
266  double /*t*/, const TimeIterator&/*interval*/) const { abort(); }
267 
271  virtual unsigned int GetControlNDoFs() const { abort(); }
279  virtual unsigned int GetStateNDoFs(int time_point = -1) const = 0;
284  virtual unsigned int GetConstraintNDoFs(std::string /*name*/) const { abort(); return 0; }
289  virtual const std::vector<unsigned int>& GetControlDoFsPerBlock() const { abort(); }
297  virtual const std::vector<unsigned int>& GetStateDoFsPerBlock(int time_point = -1) const = 0;
302  virtual const std::vector<unsigned int>& GetConstraintDoFsPerBlock(std::string /*name*/) const { abort(); }
306  virtual unsigned int GetNGlobalConstraints() const { abort(); }
310  virtual unsigned int GetNLocalConstraints() const { abort(); }
311 
315  double GetStepSize() const
316  {
317  return _interval.get_k();
318  }
322  double GetNextStepSize() const
323  {
324  assert(_interval!= _tdfh.last_interval());
325  double k = (++_interval).get_k();
326  --_interval;
327  return k;
328  }
332  double GetPreviousStepSize() const
333  {
334  assert(_interval!= _tdfh.first_interval());
335  double k = (--_interval).get_k();
336  ++_interval;
337  return k;
338  }
339 
340 
353  virtual unsigned int NewTimePointToOldTimePoint(unsigned int t) const = 0;
354 
355 
365  virtual void SpatialMeshTransferControl(const VECTOR& /*old_values*/, VECTOR& /*new_values*/) const { abort(); }
366 
367  virtual void SpatialMeshTransferState(const VECTOR& /*old_values*/, VECTOR& /*new_values*/) const { abort(); }
368 
369  protected:
374  {
375  assert( _state_ticket < std::numeric_limits<unsigned int>::max());
376  _state_ticket++;
377  }
378 
383  {
384  assert( _control_ticket < std::numeric_limits<unsigned int>::max());
385  _control_ticket++;
386  }
387 
388  private:
389  mutable TimeDoFHandler _tdfh;//FIXME Is it really necessary for _tdfh and _interval to be mutable? this is really ugly
390  mutable TimeIterator _interval;
391 
392  unsigned int _control_ticket;
393  unsigned int _state_ticket;
394  mutable DOpEtypes::ControlType _control_type;
395 };
396 
397 }
398 
399 #endif
void ReInitTime()
Definition: spacetimehandler_base.h:82
virtual unsigned int NewTimePointToOldTimePoint(unsigned int t) const =0
double GetTime(unsigned int time_point) const
Definition: spacetimehandler_base.h:139
Definition: timedofhandler.h:64
bool IsValidControlTicket(unsigned int &ticket) const
Definition: spacetimehandler_base.h:208
virtual unsigned int GetNLocalConstraints() const
Definition: spacetimehandler_base.h:310
double GetStepSize() const
Definition: spacetimehandler_base.h:315
virtual const std::vector< unsigned int > & GetConstraintDoFsPerBlock(std::string) const
Definition: spacetimehandler_base.h:302
virtual void SpatialMeshTransferState(const VECTOR &, VECTOR &) const
Definition: spacetimehandler_base.h:367
virtual void InterpolateConstraint(VECTOR &, const std::vector< VECTOR * > &, double, const TimeIterator &) const
Definition: spacetimehandler_base.h:265
virtual unsigned int GetControlNDoFs() const
Definition: spacetimehandler_base.h:271
ControlType
Definition: dopetypes.h:102
bool IsValidStateTicket(unsigned int &ticket) const
Definition: spacetimehandler_base.h:191
virtual unsigned int GetNGlobalConstraints() const
Definition: spacetimehandler_base.h:306
const std::vector< double > & GetTimes() const
Definition: timedofhandler.h:164
TimeIterator first_interval() const
Definition: timedofhandler.h:203
Definition: spacetimehandler_base.h:47
SpaceTimeHandlerBase(const dealii::Triangulation< 1 > &times, const dealii::FiniteElement< 1 > &fe, DOpEtypes::ControlType type=DOpEtypes::undefined)
Definition: spacetimehandler_base.h:64
double GetNextStepSize() const
Definition: spacetimehandler_base.h:322
SpaceTimeHandlerBase(DOpEtypes::ControlType control_type=DOpEtypes::undefined)
Definition: spacetimehandler_base.h:51
Definition: timeiterator.h:63
double get_k() const
Definition: timeiterator.h:195
unsigned int GetNbrOfDoFs() const
Definition: timedofhandler.h:152
TimeIterator last_interval() const
Definition: timedofhandler.h:224
SpaceTimeHandlerBase(const dealii::Triangulation< 1 > &times, DOpEtypes::ControlType type=DOpEtypes::undefined)
Definition: spacetimehandler_base.h:57
unsigned int GetNbrOfIntervals() const
Definition: timedofhandler.h:140
virtual void InterpolateControl(VECTOR &, const std::vector< VECTOR * > &, double, const TimeIterator &) const
Definition: spacetimehandler_base.h:238
virtual void SpatialMeshTransferControl(const VECTOR &, VECTOR &) const
Definition: spacetimehandler_base.h:365
DOpEtypes::ControlType GetControlType() const
Definition: spacetimehandler_base.h:220
const TimeDoFHandler & GetTimeDoFHandler() const
Definition: spacetimehandler_base.h:149
virtual ~SpaceTimeHandlerBase()
Definition: spacetimehandler_base.h:74
const TimeIterator & GetInterval() const
Definition: spacetimehandler_base.h:127
void SetInterval(const TimeIterator &it)
Definition: spacetimehandler_base.h:116
void distribute_dofs()
Definition: timedofhandler.h:115
virtual const std::vector< unsigned int > & GetStateDoFsPerBlock(int time_point=-1) const =0
double GetTime(unsigned int timestep)
Definition: timedofhandler.h:193
void GetTimes(const TimeIterator &interval, std::vector< double > &local_times) const
Definition: spacetimehandler_base.h:175
unsigned int GetMaxTimePoint() const
Definition: spacetimehandler_base.h:97
const std::vector< double > & GetTimes() const
Definition: spacetimehandler_base.h:158
virtual const std::vector< unsigned int > & GetControlDoFsPerBlock() const
Definition: spacetimehandler_base.h:289
void IncrementControlTicket()
Definition: spacetimehandler_base.h:382
virtual unsigned int GetStateNDoFs(int time_point=-1) const =0
double GetPreviousStepSize() const
Definition: spacetimehandler_base.h:332
virtual void InterpolateState(VECTOR &, const std::vector< VECTOR * > &, double, const TimeIterator &) const =0
unsigned int GetNbrOfIntervals() const
Definition: spacetimehandler_base.h:105
Definition: dopetypes.h:104
virtual unsigned int GetConstraintNDoFs(std::string) const
Definition: spacetimehandler_base.h:284
void IncrementStateTicket()
Definition: spacetimehandler_base.h:373