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 <deal.II/dofs/dof_handler.h>
35 #include <deal.II/dofs/dof_renumbering.h>
36 #include <deal.II/dofs/dof_tools.h>
37 #include <deal.II/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  bool flux_pattern = false,
60  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
62  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
63  dealdim>(index_setter), sparse_mkr_dynamic_(true), triangulation_(
64  triangulation), state_dof_handler_(triangulation_), state_fe_(
65  &state_fe), mapping_(
66  &DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1), state_mesh_transfer_(
67  NULL)
68  {
69  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
70  user_defined_dof_constr_ = NULL;
71  }
73  dealii::Triangulation<dealdim>& triangulation, const FE<dealdim, dealdim>& state_fe,
74  dealii::Triangulation<1> & times,
75  bool flux_pattern = false,
76  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
78  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
79  dealdim>(times, index_setter), sparse_mkr_dynamic_(true), triangulation_(
80  triangulation), state_dof_handler_(triangulation_), state_fe_(
81  &state_fe), mapping_(
82  &DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1), state_mesh_transfer_(
83  NULL)
84  {
85  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
86  user_defined_dof_constr_ = NULL;
87  }
88 
90  dealii::Triangulation<dealdim>& triangulation,
91  const DOpEWrapper::Mapping<dealdim, DH>& mapping,
92  const FE<dealdim, dealdim>& state_fe,
93  bool flux_pattern = false,
94  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
96  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
97  dealdim>(index_setter), sparse_mkr_dynamic_(true), triangulation_(
98  triangulation), state_dof_handler_(triangulation_), state_fe_(
99  &state_fe), mapping_(&mapping), state_mesh_transfer_(NULL)
100  {
101  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
102  user_defined_dof_constr_ = NULL;
103  }
105  dealii::Triangulation<dealdim>& triangulation,
106  const DOpEWrapper::Mapping<dealdim, DH>& mapping,
107  const FE<dealdim, dealdim>& state_fe,
108  dealii::Triangulation<1> & times,
109  bool flux_pattern = false,
110  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
112  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
113  dealdim>(times, index_setter), sparse_mkr_dynamic_(true), triangulation_(
114  triangulation), state_dof_handler_(triangulation_), state_fe_(
115  &state_fe), mapping_(&mapping), state_mesh_transfer_(NULL)
116  {
117  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
118  user_defined_dof_constr_ = NULL;
119  }
120 
121  virtual
123  {
124  state_dof_handler_.clear();
125 
126  if (state_mesh_transfer_ != NULL)
127  {
128  delete state_mesh_transfer_;
129  }
130  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_ == true)
131  {
132  delete sparsitymaker_;
133  }
134  }
135 
139  void
140  ReInit(unsigned int state_n_blocks,
141  const std::vector<unsigned int>& state_block_component)
142  {
143 
145  state_dof_handler_);
146  state_dof_handler_.distribute_dofs(GetFESystem("state"));
147 // DoFRenumbering::Cuthill_McKee(
148 // static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
149  DoFRenumbering::component_wise(
150  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
151 
152  state_dof_constraints_.clear();
153  DoFTools::make_hanging_node_constraints(
154  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
155  state_dof_constraints_);
156  //TODO Dirichlet ueber Constraints
157  if (GetUserDefinedDoFConstraints() != NULL
158  )
159  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
160  state_dof_handler_, state_dof_constraints_);
161  state_dof_constraints_.close();
162  state_dofs_per_block_.resize(state_n_blocks);
163 
164  DoFTools::count_dofs_per_block(
165  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
166  state_dofs_per_block_, state_block_component);
167 
168  support_points_.clear();
169 
170  //Initialize also the timediscretization.
171  this->ReInitTime();
172 
173  //There where changes invalidate tickets
174  this->IncrementStateTicket();
175  }
176 
182  {
183  //There is only one mesh, hence always return this
184  return state_dof_handler_;
185  }
186 
191  GetMapping() const
192  {
193  return *mapping_;
194  }
195 
199  unsigned int
200  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
201  {
202  return state_dofs_per_block_[b];
203  }
204 
208  const std::vector<unsigned int>&
209  GetStateDoFsPerBlock(int /*time_point*/= -1) const
210  {
211  return state_dofs_per_block_;
212  }
213 
217  const dealii::ConstraintMatrix&
219  {
220  return state_dof_constraints_;
221  }
222 
226  virtual void
227  InterpolateState(VECTOR& result,
228  const std::vector<VECTOR*> & local_vectors, double t,
229  const TimeIterator& it) const
230  {
231  assert(it.get_left() <= t);
232  assert(it.get_right() >= t);
233  if (local_vectors.size() != 2)
234  throw DOpEException(
235  "This function is currently not implemented for anything other than"
236  " linear interpolation of 2 DoFs.",
237  "MethodOfLine_SpaceTimeHandler::InterpolateState");
238 
239  double lambda_l = (it.get_right() - t) / it.get_k();
240  double lambda_r = (t - it.get_left()) / it.get_k();
241 
242  //Here we assume that the numbering of dofs goes from left to right!
243  result = *local_vectors[0];
244 
245  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
246  }
247 
251  unsigned int
252  GetStateNDoFs(int /*time_point*/= -1) const
253  {
254  return GetStateDoFHandler().n_dofs();
255  }
256 
260  const std::vector<Point<dealdim> >&
262  {
263  support_points_.resize(GetStateNDoFs());
265  < std::vector<Point<dealdim> >, dealdim
266  > (this->GetMapping(), GetStateDoFHandler(), support_points_);
267  return support_points_;
268  }
269 
270  /******************************************************/
271  void
272  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
273  {
274  this->GetSparsityMaker()->ComputeSparsityPattern(
275  this->GetStateDoFHandler(), sparsity,
277  }
278 
279  /******************************************************/
280 
284  const FE<dealdim, dealdim>&
285  GetFESystem(std::string name) const
286  {
287  if (name == "state")
288  {
289  return *state_fe_;
290  }
291  else
292  {
293  abort();
294  throw DOpEException("Not implemented for name =" + name,
295  "MethodOfLines_StateSpaceTimeHandler::GetFESystem");
296  }
297 
298  }
299 
300  /******************************************************/
311  void
314  {
315  //assert(ref_type == DOpEtypes::RefinementType::global);
316  RefinementContainer ref_con_dummy;
317  RefineSpace(ref_con_dummy);
318  }
319 
320  /******************************************************/
330  void
331  RefineSpace(const RefinementContainer& ref_container)
332  {
333  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
334 
335  //make sure that we do not use any coarsening
336  assert(!ref_container.UsesCoarsening());
337 
338  if (state_mesh_transfer_ != NULL)
339  {
340  delete state_mesh_transfer_;
341  state_mesh_transfer_ = NULL;
342  }
343  state_mesh_transfer_ = new dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >(
344  state_dof_handler_);
345 
346  if (DOpEtypes::RefinementType::global == ref_type)
347  {
348  triangulation_.set_all_refine_flags();
349  }
350  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
351  {
352  GridRefinement::refine_and_coarsen_fixed_number(triangulation_,
353  ref_container.GetLocalErrorIndicators(),
354  ref_container.GetTopFraction(),
355  ref_container.GetBottomFraction());
356  }
357  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
358  {
359 
360  GridRefinement::refine_and_coarsen_fixed_fraction(triangulation_,
361  ref_container.GetLocalErrorIndicators(),
362  ref_container.GetTopFraction(),
363  ref_container.GetBottomFraction());
364  }
365  else if (DOpEtypes::RefinementType::optimized == ref_type)
366  {
367 
368  GridRefinement::refine_and_coarsen_optimize(triangulation_,
369  ref_container.GetLocalErrorIndicators(),
370  ref_container.GetConvergenceOrder());
371  }
372  else
373  {
374  throw DOpEException("Not implemented for name =" + DOpEtypesToString(ref_type),
375  "MethodOfLines_StateSpaceTimeHandler::RefineStateSpace");
376  }
377  triangulation_.prepare_coarsening_and_refinement();
378  if (state_mesh_transfer_ != NULL
379  )
380  state_mesh_transfer_->prepare_for_pure_refinement();
381  triangulation_.execute_coarsening_and_refinement();
382  }
383  /******************************************************/
384 
388  unsigned int
389  NewTimePointToOldTimePoint(unsigned int t) const
390  {
391  //TODO this has to be implemented when temporal refinement is possible!
392  //At present the temporal grid can't be refined
393  return t;
394  }
395 
396  /******************************************************/
397 
401  void
402  SpatialMeshTransferState(const VECTOR& old_values,
403  VECTOR& new_values) const
404  {
405  if (state_mesh_transfer_ != NULL
406  )
407  state_mesh_transfer_->refine_interpolate(old_values, new_values);
408  }
409 
410  /******************************************************/
419  void
421  UserDefinedDoFConstraints<DH, dealdim>& user_defined_dof_constr)
422  {
423  user_defined_dof_constr_ = &user_defined_dof_constr;
424  user_defined_dof_constr_->RegisterMapping(this->GetMapping());
425  }
426  /******************************************************/
434  void
436  {
437  assert(sparse_mkr_dynamic_==true); //If not true, we already set the sparsity maker
438  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_)
439  delete sparsitymaker_;
440  sparsitymaker_ = &sparsity_maker;
441  sparse_mkr_dynamic_ = false;
442  }
443 
444  private:
446  GetSparsityMaker() const
447  {
448  return sparsitymaker_;
449  }
450  const UserDefinedDoFConstraints<DH, dealdim>*
451  GetUserDefinedDoFConstraints() const
452  {
453  return user_defined_dof_constr_;
454  }
455  SparsityMaker<DH, dealdim>* sparsitymaker_;
456  UserDefinedDoFConstraints<DH, dealdim>* user_defined_dof_constr_;
457  bool sparse_mkr_dynamic_;
458 
459  dealii::Triangulation<dealdim>& triangulation_;
460  DOpEWrapper::DoFHandler<dealdim, DH> state_dof_handler_;
461 
462  std::vector<unsigned int> state_dofs_per_block_;
463 
464  dealii::ConstraintMatrix state_dof_constraints_;
465 
466  const dealii::SmartPointer<const FE<dealdim, dealdim> > state_fe_; //TODO is there a reason that this is not a reference?
467  const dealii::SmartPointer<
468  const DOpEWrapper::Mapping<dealdim, DH> > mapping_;
469 
470  std::vector<Point<dealdim> > support_points_;
471  dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >* state_mesh_transfer_;
472 
473  };
474 
475 }
476 #endif
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_statespacetimehandler.h:218
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:140
void RefineSpace(DOpEtypes::RefinementType=DOpEtypes::RefinementType::global)
Definition: mol_statespacetimehandler.h:312
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &state_fe, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:57
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: statespacetimehandler.h:253
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:72
Definition: dopetypes.h:52
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_statespacetimehandler.h:389
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_statespacetimehandler.h:331
Definition: mol_statespacetimehandler.h:53
const DOpEWrapper::Mapping< dealdim, DH > & GetMapping() const
Definition: mol_statespacetimehandler.h:191
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_statespacetimehandler.h:252
Definition: timeiterator.h:62
double get_k() const
Definition: timeiterator.h:194
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:86
RefinementType
Definition: dopetypes.h:50
void SetSparsityMaker(SparsityMaker< DH, dealdim > &sparsity_maker)
Definition: mol_statespacetimehandler.h:435
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_statespacetimehandler.h:272
Definition: refinementcontainer.h:42
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:122
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:76
const DOpEWrapper::DoFHandler< dealdim, DH > & GetStateDoFHandler() const
Definition: mol_statespacetimehandler.h:181
double get_right() const
Definition: timeiterator.h:184
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:420
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_statespacetimehandler.h:402
const FE< dealdim, dealdim > & GetFESystem(std::string name) const
Definition: mol_statespacetimehandler.h:285
const std::vector< Point< dealdim > > & GetMapDoFToSupportPoints()
Definition: mol_statespacetimehandler.h:261
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:113
unsigned int GetStateDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_statespacetimehandler.h:200
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const DOpEWrapper::Mapping< dealdim, DH > &mapping, const FE< dealdim, dealdim > &state_fe, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:89
bool UsesCoarsening() const
Definition: refinementcontainer.cc:94
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_statespacetimehandler.h:209
double get_left() const
Definition: timeiterator.h:174
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const DOpEWrapper::Mapping< dealdim, DH > &mapping, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:104
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:227
Definition: dopetypes.h:52
virtual double GetBottomFraction() const
Definition: refinementcontainer.cc:66
std::string DOpEtypesToString(const C &)
Definition: dopetypes.h:134
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