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  const 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, const dealii::Triangulation<1> & times,
105  const ActiveFEIndexSetterInterface<dealdim>& index_setter =
107  : StateSpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
108  dealdim>(times, index_setter), _sparse_mkr_dynamic(true), _triangulation(
109  triangulation), _state_dof_handler(_triangulation), _state_fe(
110  &state_fe), _mapping(&mapping), _state_mesh_transfer(NULL)
111  {
112  _sparsitymaker = new SparsityMaker<DH, dealdim>;
113  _user_defined_dof_constr = NULL;
114  }
115 
116  virtual
118  {
119  _state_dof_handler.clear();
120 
121  if (_state_mesh_transfer != NULL)
122  {
123  delete _state_mesh_transfer;
124  }
125  if (_sparsitymaker != NULL && _sparse_mkr_dynamic == true)
126  {
127  delete _sparsitymaker;
128  }
129  }
130 
134  void
135  ReInit(unsigned int state_n_blocks,
136  const std::vector<unsigned int>& state_block_component)
137  {
138 
140  _state_dof_handler);
141  _state_dof_handler.distribute_dofs(GetFESystem("state"));
142  DoFRenumbering::Cuthill_McKee(
143  static_cast<DH<dealdim, dealdim>&>(_state_dof_handler));
144  DoFRenumbering::component_wise(
145  static_cast<DH<dealdim, dealdim>&>(_state_dof_handler));
146 
147  _state_dof_constraints.clear();
148  DoFTools::make_hanging_node_constraints(
149  static_cast<DH<dealdim, dealdim>&>(_state_dof_handler),
150  _state_dof_constraints);
151  //TODO Dirichlet ueber Constraints
152  if (GetUserDefinedDoFConstraints() != NULL
153  )
154  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
155  _state_dof_handler, _state_dof_constraints);
156  _state_dof_constraints.close();
157  _state_dofs_per_block.resize(state_n_blocks);
158 
159  DoFTools::count_dofs_per_block(
160  static_cast<DH<dealdim, dealdim>&>(_state_dof_handler),
161  _state_dofs_per_block, state_block_component);
162 
163  _support_points.clear();
164 
165  //Initialize also the timediscretization.
166  this->ReInitTime();
167 
168  //There where changes invalidate tickets
169  this->IncrementStateTicket();
170  }
171 
177  {
178  //There is only one mesh, hence always return this
179  return _state_dof_handler;
180  }
181 
186  GetMapping() const
187  {
188  return *_mapping;
189  }
190 
194  unsigned int
195  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
196  {
197  return _state_dofs_per_block[b];
198  }
199 
203  const std::vector<unsigned int>&
204  GetStateDoFsPerBlock(int /*time_point*/= -1) const
205  {
206  return _state_dofs_per_block;
207  }
208 
212  const dealii::ConstraintMatrix&
214  {
215  return _state_dof_constraints;
216  }
217 
221  virtual void
222  InterpolateState(VECTOR& result,
223  const std::vector<VECTOR*> & local_vectors, double t,
224  const TimeIterator& it) const
225  {
226  assert(it.get_left() <= t);
227  assert(it.get_right() >= t);
228  if (local_vectors.size() != 2)
229  throw DOpEException(
230  "This function is currently not implemented for anything other than"
231  " linear interpolation of 2 DoFs.",
232  "MethodOfLine_SpaceTimeHandler::InterpolateState");
233 
234  double lambda_l = (it.get_right() - t) / it.get_k();
235  double lambda_r = (t - it.get_left()) / it.get_k();
236 
237  //Here we assume that the numbering of dofs goes from left to right!
238  result = *local_vectors[0];
239 
240  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
241  }
242 
246  unsigned int
247  GetStateNDoFs(int /*time_point*/= -1) const
248  {
249  return GetStateDoFHandler().n_dofs();
250  }
251 
255  const std::vector<Point<dealdim> >&
257  {
258  _support_points.resize(GetStateNDoFs());
260  < std::vector<Point<dealdim> >, dealdim
261  > (this->GetMapping(), GetStateDoFHandler(), _support_points);
262  return _support_points;
263  }
264 
265  /******************************************************/
266  void
267  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
268  {
269  this->GetSparsityMaker()->ComputeSparsityPattern(
270  this->GetStateDoFHandler(), sparsity,
272  }
273 
274  /******************************************************/
275 
279  const FE<dealdim, dealdim>&
280  GetFESystem(std::string name) const
281  {
282  if (name == "state")
283  {
284  return *_state_fe;
285  }
286  else
287  {
288  abort();
289  throw DOpEException("Not implemented for name =" + name,
290  "MethodOfLines_StateSpaceTimeHandler::GetFESystem");
291  }
292 
293  }
294 
295  /******************************************************/
306  void
309  {
310  assert(ref_type == DOpEtypes::RefinementType::global);
311  RefinementContainer ref_con_dummy;
312  RefineSpace(ref_con_dummy);
313  }
314 
315  /******************************************************/
325  void
326  RefineSpace(const RefinementContainer& ref_container)
327  {
328  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
329 
330  //make sure that we do not use any coarsening
331  assert(!ref_container.UsesCoarsening());
332 
333  if (_state_mesh_transfer != NULL)
334  {
335  delete _state_mesh_transfer;
336  _state_mesh_transfer = NULL;
337  }
338  _state_mesh_transfer = new dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >(
339  _state_dof_handler);
340 
341  if (DOpEtypes::RefinementType::global == ref_type)
342  {
343  _triangulation.set_all_refine_flags();
344  }
345  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
346  {
347  GridRefinement::refine_and_coarsen_fixed_number(_triangulation,
348  ref_container.GetLocalErrorIndicators(),
349  ref_container.GetTopFraction(),
350  ref_container.GetBottomFraction());
351  }
352  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
353  {
354 
355  GridRefinement::refine_and_coarsen_fixed_fraction(_triangulation,
356  ref_container.GetLocalErrorIndicators(),
357  ref_container.GetTopFraction(),
358  ref_container.GetBottomFraction());
359  }
360  else if (DOpEtypes::RefinementType::optimized == ref_type)
361  {
362 
363  GridRefinement::refine_and_coarsen_optimize(_triangulation,
364  ref_container.GetLocalErrorIndicators(),
365  ref_container.GetConvergenceOrder());
366  }
367  else
368  {
369  throw DOpEException("Not implemented for name =" + ref_type,
370  "MethodOfLines_StateSpaceTimeHandler::RefineStateSpace");
371  }
372  _triangulation.prepare_coarsening_and_refinement();
373  if (_state_mesh_transfer != NULL
374  )
375  _state_mesh_transfer->prepare_for_pure_refinement();
376  _triangulation.execute_coarsening_and_refinement();
377  }
378  /******************************************************/
379 
383  unsigned int
384  NewTimePointToOldTimePoint(unsigned int t) const
385  {
386  //TODO this has to be implemented when temporal refinement is possible!
387  //At present the temporal grid can't be refined
388  return t;
389  }
390 
391  /******************************************************/
392 
396  void
397  SpatialMeshTransferState(const VECTOR& old_values,
398  VECTOR& new_values) const
399  {
400  if (_state_mesh_transfer != NULL
401  )
402  _state_mesh_transfer->refine_interpolate(old_values, new_values);
403  }
404 
405  /******************************************************/
414  void
416  UserDefinedDoFConstraints<DH, dealdim>& user_defined_dof_constr)
417  {
418  _user_defined_dof_constr = &user_defined_dof_constr;
419  _user_defined_dof_constr->RegisterMapping(this->GetMapping());
420  }
421  /******************************************************/
429  void
431  {
432  if (_sparsitymaker != NULL && _sparse_mkr_dynamic)
433  delete _sparsitymaker;
434  _sparsitymaker = &sparsity_maker;
435  _sparse_mkr_dynamic = false;
436  }
437 
438  private:
440  GetSparsityMaker() const
441  {
442  return _sparsitymaker;
443  }
444  const UserDefinedDoFConstraints<DH, dealdim>*
445  GetUserDefinedDoFConstraints() const
446  {
447  return _user_defined_dof_constr;
448  }
449  SparsityMaker<DH, dealdim>* _sparsitymaker;
450  UserDefinedDoFConstraints<DH, dealdim>* _user_defined_dof_constr;
451  bool _sparse_mkr_dynamic;
452 
453  dealii::Triangulation<dealdim>& _triangulation;
454  DOpEWrapper::DoFHandler<dealdim, DH> _state_dof_handler;
455 
456  std::vector<unsigned int> _state_dofs_per_block;
457 
458  dealii::ConstraintMatrix _state_dof_constraints;
459 
460  const dealii::SmartPointer<const FE<dealdim, dealdim> > _state_fe; //TODO is there a reason that this is not a reference?
461  const dealii::SmartPointer<
462  const DOpEWrapper::Mapping<dealdim, DH> > _mapping;
463 
464  std::vector<Point<dealdim> > _support_points;
465  dealii::SolutionTransfer<dealdim, VECTOR, DH<dealdim, dealdim> >* _state_mesh_transfer;
466 
467  };
468 
469 }
470 #endif
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_statespacetimehandler.h:213
void ReInitTime()
Definition: spacetimehandler_base.h:82
void RegisterMapping(const typename DOpEWrapper::Mapping< dealdim, DH > &mapping)
Definition: userdefineddofconstraints.h:76
Definition: dopetypes.h:50
void ReInit(unsigned int state_n_blocks, const std::vector< unsigned int > &state_block_component)
Definition: mol_statespacetimehandler.h:135
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:307
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: statespacetimehandler.h:253
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const DOpEWrapper::Mapping< dealdim, DH > &mapping, const FE< dealdim, dealdim > &state_fe, const dealii::Triangulation< 1 > &times, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:101
Definition: dopetypes.h:50
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_statespacetimehandler.h:384
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_statespacetimehandler.h:326
Definition: mol_statespacetimehandler.h:53
MethodOfLines_StateSpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &state_fe, const dealii::Triangulation< 1 > &times, const ActiveFEIndexSetterInterface< dealdim > &index_setter=ActiveFEIndexSetterInterface< dealdim >())
Definition: mol_statespacetimehandler.h:71
const DOpEWrapper::Mapping< dealdim, DH > & GetMapping() const
Definition: mol_statespacetimehandler.h:186
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_statespacetimehandler.h:247
Definition: timeiterator.h:63
double get_k() const
Definition: timeiterator.h:195
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:95
RefinementType
Definition: dopetypes.h:48
void SetSparsityMaker(SparsityMaker< DH, dealdim > &sparsity_maker)
Definition: mol_statespacetimehandler.h:430
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_statespacetimehandler.h:267
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:50
virtual ~MethodOfLines_StateSpaceTimeHandler()
Definition: mol_statespacetimehandler.h:117
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:85
const DOpEWrapper::DoFHandler< dealdim, DH > & GetStateDoFHandler() const
Definition: mol_statespacetimehandler.h:176
double get_right() const
Definition: timeiterator.h:185
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:415
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_statespacetimehandler.h:397
const FE< dealdim, dealdim > & GetFESystem(std::string name) const
Definition: mol_statespacetimehandler.h:280
const std::vector< Point< dealdim > > & GetMapDoFToSupportPoints()
Definition: mol_statespacetimehandler.h:256
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:195
bool UsesCoarsening() const
Definition: refinementcontainer.cc:103
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_statespacetimehandler.h:204
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:222
Definition: dopetypes.h:50
virtual double GetBottomFraction() const
Definition: refinementcontainer.cc:66
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:373