DOpE
mol_statespacetimehandler.h
Go to the documentation of this file.
1 
24 #ifndef MOL_STATESPACE_TIME_HANDLER_H_
25 #define MOL_STATESPACE_TIME_HANDLER_H_
26 
27 #include "statespacetimehandler.h"
28 #include "constraints.h"
29 #include "sparsitymaker.h"
31 #include "sth_internals.h"
32 #include "refinementcontainer.h"
33 
34 #include <dofs/dof_handler.h>
35 #include <dofs/dof_renumbering.h>
36 #include <dofs/dof_tools.h>
37 #include <lac/constraint_matrix.h>
38 #include <deal.II/numerics/solution_transfer.h>
39 #include <deal.II/grid/grid_refinement.h>
40 
41 namespace DOpE
42 {
51  template<template<int, int> class FE, template<int, int> class DH, typename SPARSITYPATTERN,
52  typename VECTOR, int dealdim>
54  DH, SPARSITYPATTERN, VECTOR, dealdim>
55  {
56  public:
58  dealii::Triangulation<dealdim>& triangulation, const FE<dealdim, dealdim>& state_fe,
59  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
61  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
62  dealdim>(index_setter), sparse_mkr_dynamic_(true), triangulation_(
63  triangulation), state_dof_handler_(triangulation_), state_fe_(
64  &state_fe), mapping_(
65  &DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1), state_mesh_transfer_(
66  NULL)
67  {
68  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
69  user_defined_dof_constr_ = NULL;
70  }
72  dealii::Triangulation<dealdim>& triangulation, const FE<dealdim, dealdim>& state_fe,
73  dealii::Triangulation<1> & times,
74  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
76  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
77  dealdim>(times, index_setter), sparse_mkr_dynamic_(true), triangulation_(
78  triangulation), state_dof_handler_(triangulation_), state_fe_(
79  &state_fe), mapping_(
80  &DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1), state_mesh_transfer_(
81  NULL)
82  {
83  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
84  user_defined_dof_constr_ = NULL;
85  }
86 
88  dealii::Triangulation<dealdim>& triangulation,
89  const DOpEWrapper::Mapping<dealdim, DH>& mapping,
90  const FE<dealdim, dealdim>& state_fe,
91  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
93  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
94  dealdim>(index_setter), sparse_mkr_dynamic_(true), triangulation_(
95  triangulation), state_dof_handler_(triangulation_), state_fe_(
96  &state_fe), mapping_(&mapping), state_mesh_transfer_(NULL)
97  {
98  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
99  user_defined_dof_constr_ = NULL;
100  }
102  dealii::Triangulation<dealdim>& triangulation,
103  const DOpEWrapper::Mapping<dealdim, DH>& mapping,
104  const FE<dealdim, dealdim>& state_fe,
105  dealii::Triangulation<1> & times,
106  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
108  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
109  dealdim>(times, index_setter), sparse_mkr_dynamic_(true), triangulation_(
110  triangulation), state_dof_handler_(triangulation_), state_fe_(
111  &state_fe), mapping_(&mapping), state_mesh_transfer_(NULL)
112  {
113  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
114  user_defined_dof_constr_ = NULL;
115  }
116 
117  virtual
119  {
120  state_dof_handler_.clear();
121 
122  if (state_mesh_transfer_ != NULL)
123  {
124  delete state_mesh_transfer_;
125  }
126  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_ == true)
127  {
128  delete sparsitymaker_;
129  }
130  }
131 
135  void
136  ReInit(unsigned int state_n_blocks,
137  const std::vector<unsigned int>& state_block_component)
138  {
139 
141  state_dof_handler_);
142  state_dof_handler_.distribute_dofs(GetFESystem("state"));
143  DoFRenumbering::Cuthill_McKee(
144  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
145  DoFRenumbering::component_wise(
146  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
147 
148  state_dof_constraints_.clear();
149  DoFTools::make_hanging_node_constraints(
150  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
151  state_dof_constraints_);
152  //TODO Dirichlet ueber Constraints
153  if (GetUserDefinedDoFConstraints() != NULL
154  )
155  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
156  state_dof_handler_, state_dof_constraints_);
157  state_dof_constraints_.close();
158  state_dofs_per_block_.resize(state_n_blocks);
159 
160  DoFTools::count_dofs_per_block(
161  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
162  state_dofs_per_block_, state_block_component);
163 
164  support_points_.clear();
165 
166  //Initialize also the timediscretization.
167  this->ReInitTime();
168 
169  //There where changes invalidate tickets
170  this->IncrementStateTicket();
171  }
172 
178  {
179  //There is only one mesh, hence always return this
180  return state_dof_handler_;
181  }
182 
187  GetMapping() const
188  {
189  return *mapping_;
190  }
191 
195  unsigned int
196  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
197  {
198  return state_dofs_per_block_[b];
199  }
200 
204  const std::vector<unsigned int>&
205  GetStateDoFsPerBlock(int /*time_point*/= -1) const
206  {
207  return state_dofs_per_block_;
208  }
209 
213  const dealii::ConstraintMatrix&
215  {
216  return state_dof_constraints_;
217  }
218 
222  virtual void
223  InterpolateState(VECTOR& result,
224  const std::vector<VECTOR*> & local_vectors, double t,
225  const TimeIterator& it) const
226  {
227  assert(it.get_left() <= t);
228  assert(it.get_right() >= t);
229  if (local_vectors.size() != 2)
230  throw DOpEException(
231  "This function is currently not implemented for anything other than"
232  " linear interpolation of 2 DoFs.",
233  "MethodOfLine_SpaceTimeHandler::InterpolateState");
234 
235  double lambda_l = (it.get_right() - t) / it.get_k();
236  double lambda_r = (t - it.get_left()) / it.get_k();
237 
238  //Here we assume that the numbering of dofs goes from left to right!
239  result = *local_vectors[0];
240 
241  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
242  }
243 
247  unsigned int
248  GetStateNDoFs(int /*time_point*/= -1) const
249  {
250  return GetStateDoFHandler().n_dofs();
251  }
252 
256  const std::vector<Point<dealdim> >&
258  {
259  support_points_.resize(GetStateNDoFs());
261  < std::vector<Point<dealdim> >, dealdim
262  > (this->GetMapping(), GetStateDoFHandler(), support_points_);
263  return support_points_;
264  }
265 
266  /******************************************************/
267  void
268  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
269  {
270  this->GetSparsityMaker()->ComputeSparsityPattern(
271  this->GetStateDoFHandler(), sparsity,
273  }
274 
275  /******************************************************/
276 
280  const FE<dealdim, dealdim>&
281  GetFESystem(std::string name) const
282  {
283  if (name == "state")
284  {
285  return *state_fe_;
286  }
287  else
288  {
289  abort();
290  throw DOpEException("Not implemented for name =" + name,
291  "MethodOfLines_StateSpaceTimeHandler::GetFESystem");
292  }
293 
294  }
295 
296  /******************************************************/
307  void
310  {
311  assert(ref_type == DOpEtypes::RefinementType::global);
312  RefinementContainer ref_con_dummy;
313  RefineSpace(ref_con_dummy);
314  }
315 
316  /******************************************************/
326  void
327  RefineSpace(const RefinementContainer& ref_container)
328  {
329  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
330 
331  //make sure that we do not use any coarsening
332  assert(!ref_container.UsesCoarsening());
333 
334  if (state_mesh_transfer_ != NULL)
335  {
336  delete state_mesh_transfer_;
337  state_mesh_transfer_ = NULL;
338  }
339  state_mesh_transfer_ = new dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >(
340  state_dof_handler_);
341 
342  if (DOpEtypes::RefinementType::global == ref_type)
343  {
344  triangulation_.set_all_refine_flags();
345  }
346  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
347  {
348  GridRefinement::refine_and_coarsen_fixed_number(triangulation_,
349  ref_container.GetLocalErrorIndicators(),
350  ref_container.GetTopFraction(),
351  ref_container.GetBottomFraction());
352  }
353  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
354  {
355 
356  GridRefinement::refine_and_coarsen_fixed_fraction(triangulation_,
357  ref_container.GetLocalErrorIndicators(),
358  ref_container.GetTopFraction(),
359  ref_container.GetBottomFraction());
360  }
361  else if (DOpEtypes::RefinementType::optimized == ref_type)
362  {
363 
364  GridRefinement::refine_and_coarsen_optimize(triangulation_,
365  ref_container.GetLocalErrorIndicators(),
366  ref_container.GetConvergenceOrder());
367  }
368  else
369  {
370  throw DOpEException("Not implemented for name =" + ref_type,
371  "MethodOfLines_StateSpaceTimeHandler::RefineStateSpace");
372  }
373  triangulation_.prepare_coarsening_and_refinement();
374  if (state_mesh_transfer_ != NULL
375  )
376  state_mesh_transfer_->prepare_for_pure_refinement();
377  triangulation_.execute_coarsening_and_refinement();
378  }
379  /******************************************************/
380 
384  unsigned int
385  NewTimePointToOldTimePoint(unsigned int t) const
386  {
387  //TODO this has to be implemented when temporal refinement is possible!
388  //At present the temporal grid can't be refined
389  return t;
390  }
391 
392  /******************************************************/
393 
397  void
398  SpatialMeshTransferState(const VECTOR& old_values,
399  VECTOR& new_values) const
400  {
401  if (state_mesh_transfer_ != NULL
402  )
403  state_mesh_transfer_->refine_interpolate(old_values, new_values);
404  }
405 
406  /******************************************************/
415  void
417  UserDefinedDoFConstraints<DH, dealdim>& user_defined_dof_constr)
418  {
419  user_defined_dof_constr_ = &user_defined_dof_constr;
420  user_defined_dof_constr_->RegisterMapping(this->GetMapping());
421  }
422  /******************************************************/
430  void
432  {
433  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_)
434  delete sparsitymaker_;
435  sparsitymaker_ = &sparsity_maker;
436  sparse_mkr_dynamic_ = false;
437  }
438 
439  private:
441  GetSparsityMaker() const
442  {
443  return sparsitymaker_;
444  }
445  const UserDefinedDoFConstraints<DH, dealdim>*
446  GetUserDefinedDoFConstraints() const
447  {
448  return user_defined_dof_constr_;
449  }
450  SparsityMaker<DH, dealdim>* sparsitymaker_;
451  UserDefinedDoFConstraints<DH, dealdim>* user_defined_dof_constr_;
452  bool sparse_mkr_dynamic_;
453 
454  dealii::Triangulation<dealdim>& triangulation_;
455  DOpEWrapper::DoFHandler<dealdim, DH> state_dof_handler_;
456 
457  std::vector<unsigned int> state_dofs_per_block_;
458 
459  dealii::ConstraintMatrix state_dof_constraints_;
460 
461  const dealii::SmartPointer<const FE<dealdim, dealdim> > state_fe_; //TODO is there a reason that this is not a reference?
462  const dealii::SmartPointer<
463  const DOpEWrapper::Mapping<dealdim, DH> > mapping_;
464 
465  std::vector<Point<dealdim> > support_points_;
466  dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >* state_mesh_transfer_;
467 
468  };
469 
470 }
471 #endif
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_statespacetimehandler.h:214
void ReInitTime()
Definition: spacetimehandler_base.h:86
void RegisterMapping(const typename DOpEWrapper::Mapping< dealdim, DH > &mapping)
Definition: userdefineddofconstraints.h:76
Definition: dopetypes.h:52
void ReInit(unsigned int state_n_blocks, const std::vector< unsigned int > &state_block_component)
Definition: mol_statespacetimehandler.h:136
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const DOpEWrapper::Mapping< dealdim, DH > &mapping, const FE< dealdim, dealdim > &state_fe, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:87
void RefineSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_statespacetimehandler.h:308
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: statespacetimehandler.h:253
Definition: dopetypes.h:52
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_statespacetimehandler.h:385
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_statespacetimehandler.h:327
Definition: mol_statespacetimehandler.h:53
const DOpEWrapper::Mapping< dealdim, DH > & GetMapping() const
Definition: mol_statespacetimehandler.h:187
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_statespacetimehandler.h:248
Definition: timeiterator.h:63
double get_k() const
Definition: timeiterator.h:195
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:95
RefinementType
Definition: dopetypes.h:50
void SetSparsityMaker(SparsityMaker< DH, dealdim > &sparsity_maker)
Definition: mol_statespacetimehandler.h:431
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_statespacetimehandler.h:268
Definition: refinementcontainer.h:42
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &state_fe, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:57
virtual const dealii::Vector< float > & GetLocalErrorIndicators() const
Definition: refinementcontainer.cc:46
virtual double GetTopFraction() const
Definition: refinementcontainer.cc:56
Definition: statespacetimehandler.h:59
Definition: dopetypes.h:52
virtual ~MethodOfLines_StateSpaceTimeHandler()
Definition: mol_statespacetimehandler.h:118
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:85
const DOpEWrapper::DoFHandler< dealdim, DH > & GetStateDoFHandler() const
Definition: mol_statespacetimehandler.h:177
double get_right() const
Definition: timeiterator.h:185
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const DOpEWrapper::Mapping< dealdim, DH > &mapping, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:101
void MapDoFsToSupportPoints(const DOpEWrapper::Mapping< dealdim, dealii::DoFHandler > &mapping, const DOpEWrapper::DoFHandler< dealdim, dealii::DoFHandler > &dh, VECTOR &support_points)
Definition: sth_internals.h:47
void SetUserDefinedDoFConstraints(UserDefinedDoFConstraints< DH, dealdim > &user_defined_dof_constr)
Definition: mol_statespacetimehandler.h:416
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_statespacetimehandler.h:398
const FE< dealdim, dealdim > & GetFESystem(std::string name) const
Definition: mol_statespacetimehandler.h:281
const std::vector< Point< dealdim > > & GetMapDoFToSupportPoints()
Definition: mol_statespacetimehandler.h:257
virtual void ComputeSparsityPattern(const DOpEWrapper::DoFHandler< dim, DH > &dof_handler, dealii::BlockSparsityPattern &sparsity, const dealii::ConstraintMatrix &hanging_node_constraints, const std::vector< unsigned int > &blocks) const
Definition: sparsitymaker.h:111
unsigned int GetStateDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_statespacetimehandler.h:196
bool UsesCoarsening() const
Definition: refinementcontainer.cc:103
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_statespacetimehandler.h:205
double get_left() const
Definition: timeiterator.h:175
Definition: dopeexception.h:35
virtual void InterpolateState(VECTOR &result, const std::vector< VECTOR * > &local_vectors, double t, const TimeIterator &it) const
Definition: mol_statespacetimehandler.h:223
Definition: dopetypes.h:52
virtual double GetBottomFraction() const
Definition: refinementcontainer.cc:66
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:71
virtual void MakeStateDoFConstraints(const DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler, dealii::ConstraintMatrix &dof_constraints) const
Definition: userdefineddofconstraints.h:93
void IncrementStateTicket()
Definition: spacetimehandler_base.h:437