DOpE
mol_spacetimehandler.h
Go to the documentation of this file.
1 
24 #ifndef MOL_SPACE_TIME_HANDLER_H_
25 #define MOL_SPACE_TIME_HANDLER_H_
26 
27 #include "spacetimehandler.h"
28 #include "constraints.h"
29 #include "sparsitymaker.h"
31 #include "sth_internals.h"
32 #include "mapping_wrapper.h"
33 #include "refinementcontainer.h"
35 
36 #include <dofs/dof_handler.h>
37 #include <dofs/dof_renumbering.h>
38 #include <dofs/dof_tools.h>
39 #include <hp/mapping_collection.h>
40 #include <lac/constraint_matrix.h>
41 #include <deal.II/grid/grid_refinement.h>
42 
43 namespace DOpE
44 {
49  template<template <int, int> class FE, template<int, int> class DH, typename SPARSITYPATTERN,
50  typename VECTOR, int dopedim, int dealdim>
52  SPARSITYPATTERN, VECTOR, dopedim, dealdim>
53  {
54  public:
65  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
66  const FE<dealdim, dealdim>& control_fe,
67  const FE<dealdim, dealdim>& state_fe,
71  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(type, index_setter),
72  triangulation_(triangulation),
73  control_dof_handler_(triangulation_),
74  state_dof_handler_(triangulation_),
75  control_fe_(&control_fe),
76  state_fe_(&state_fe),
77  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
78  constraints_(),
79  control_mesh_transfer_(NULL),
80  state_mesh_transfer_(NULL),
81  sparse_mkr_dynamic_(true)
82  {
83  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
84  user_defined_dof_constr_ = NULL;
85  }
86 
98  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
99  const FE<dealdim, dealdim>& control_fe,
100  const FE<dealdim, dealdim>& state_fe,
101  dealii::Triangulation<1> & times,
105  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(times, type, index_setter),
106  triangulation_(triangulation),
107  control_dof_handler_(triangulation_),
108  state_dof_handler_(triangulation_),
109  control_fe_(&control_fe),
110  state_fe_(&state_fe),
111  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
112  constraints_(),
113  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
114  {
115  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
116  user_defined_dof_constr_ = NULL;
117  }
118 
130  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
131  const FE<dealdim, dealdim>& control_fe,
132  const FE<dealdim, dealdim>& state_fe,
133  const Constraints& c,
137  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(type, index_setter),
138  triangulation_(triangulation),
139  control_dof_handler_(triangulation_),
140  state_dof_handler_(triangulation_),
141  control_fe_(&control_fe),
142  state_fe_(&state_fe),
143  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
144  constraints_(c),
145  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
146  {
147  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
148  user_defined_dof_constr_ = NULL;
149  }
150 
163  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
164  const FE<dealdim, dealdim>& control_fe,
165  const FE<dealdim, dealdim>& state_fe,
166  dealii::Triangulation<1> & times,
167  const Constraints& c,
171  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(times, type, index_setter),
172  triangulation_(triangulation),
173  control_dof_handler_(triangulation_),
174  state_dof_handler_(triangulation_),
175  control_fe_(&control_fe),
176  state_fe_(&state_fe),
177  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
178  constraints_(c),
179  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
180  {
181  sparsitymaker_ = new SparsityMaker<DH, dealdim>;
182  user_defined_dof_constr_ = NULL;
183  }
184 
186  {
187  control_dof_handler_.clear();
188 
189  state_dof_handler_.clear();
190 
191  if (control_mesh_transfer_ != NULL)
192  {
193  delete control_mesh_transfer_;
194  }
195  if (state_mesh_transfer_ != NULL)
196  {
197  delete state_mesh_transfer_;
198  }
199  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_ == true)
200  {
201  delete sparsitymaker_;
202  }
203  }
204 
208  void
209  ReInit(unsigned int control_n_blocks,
210  const std::vector<unsigned int>& control_block_component,
211  unsigned int state_n_blocks,
212  const std::vector<unsigned int>& state_block_component)
213  {
214 #if dope_dimension > 0
215  SpaceTimeHandler<FE, DH, SPARSITYPATTERN,
216  VECTOR, dopedim, dealdim>::SetActiveFEIndicesControl(control_dof_handler_);
217 #endif
218  control_dof_handler_.distribute_dofs(*control_fe_);
219 
220 #if dope_dimension > 0
221  DoFRenumbering::component_wise (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_));
222  if(dopedim==dealdim)
223  {
224  control_dof_constraints_.clear ();
225  DoFTools::make_hanging_node_constraints (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_),
226  control_dof_constraints_);
227  if (GetUserDefinedDoFConstraints() != NULL)
228  GetUserDefinedDoFConstraints()->MakeControlDoFConstraints(control_dof_handler_,
229  control_dof_constraints_);
230  control_dof_constraints_.close ();
231  }
232  else
233  {
234  throw DOpEException("Not implemented for dopedim != dealdim","MethodOfLines_SpaceTimeHandler::ReInit");
235  }
236 #endif
237  control_dofs_per_block_.resize(control_n_blocks);
238 #if dope_dimension > 0
239  {
240  DoFTools::count_dofs_per_block (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_),
241  control_dofs_per_block_,control_block_component);
242  }
243 #else
244  {
245  for (unsigned int i = 0; i < control_dofs_per_block_.size(); i++)
246  {
247  control_dofs_per_block_[i] = 0;
248  }
249  for (unsigned int i = 0; i < control_block_component.size(); i++)
250  {
251  control_dofs_per_block_[control_block_component[i]]++;
252  }
253  }
254 #endif
256  state_dof_handler_);
257  state_dof_handler_.distribute_dofs(GetFESystem("state"));
258  DoFRenumbering::component_wise(
259  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
260 
261  state_dof_constraints_.clear();
262  DoFTools::make_hanging_node_constraints(
263  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
264  state_dof_constraints_);
265  //TODO Dirichlet ueber Constraints
266  if (GetUserDefinedDoFConstraints() != NULL)
267  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
268  state_dof_handler_, state_dof_constraints_);
269  state_dof_constraints_.close();
270 
271  state_dofs_per_block_.resize(state_n_blocks);
272  DoFTools::count_dofs_per_block(
273  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
274  state_dofs_per_block_, state_block_component);
275 
276  support_points_.clear();
277 
278  constraints_.ReInit(control_dofs_per_block_);
279  //constraints_.ReInit(control_dofs_per_block_, state_dofs_per_block_);
280 
281  //Initialize also the timediscretization.
282  this->ReInitTime();
283 
284  //There where changes invalidate tickets
285  this->IncrementControlTicket();
286  this->IncrementStateTicket();
287  }
288 
294  {
295  //There is only one mesh, hence always return this
296  return control_dof_handler_;
297  }
303  {
304  //There is only one mesh, hence always return this
305  return state_dof_handler_;
306  }
311  GetMapping() const
312  {
313  return *mapping_;
314  }
315 
319  unsigned int
320  GetControlDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
321  {
322  return control_dofs_per_block_[b];
323  }
327  unsigned int
328  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
329  {
330  return state_dofs_per_block_[b];
331  }
335  unsigned int
336  GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
337  {
338  return (constraints_.GetDoFsPerBlock(name))[b];
339  }
343  const std::vector<unsigned int>&
344  GetControlDoFsPerBlock(int /*time_point*/= -1) const
345  {
346  return control_dofs_per_block_;
347  }
351  const std::vector<unsigned int>&
352  GetStateDoFsPerBlock(int /*time_point*/= -1) const
353  {
354  return state_dofs_per_block_;
355  }
359  const std::vector<unsigned int>&
360  GetConstraintDoFsPerBlock(std::string name) const
361  {
362  return constraints_.GetDoFsPerBlock(name);
363  }
367  const dealii::ConstraintMatrix&
369  {
370  return control_dof_constraints_;
371  }
375  const dealii::ConstraintMatrix&
377  {
378  return state_dof_constraints_;
379  }
380 
384  virtual void
385  InterpolateControl(VECTOR& result,
386  const std::vector<VECTOR*> & local_vectors, double t,
387  const TimeIterator& it) const
388  {
389  assert(it.get_left() <= t);
390  assert(it.get_right() >= t);
391  if (local_vectors.size() != 2)
392  {
393  throw DOpEException(
394  "This function is currently not implemented for anything other than"
395  " linear interpolation of 2 DoFs.",
396  "MethodOfLine_SpaceTimeHandler::InterpolateControl");
397  }
398  double lambda_l = (it.get_right() - t) / it.get_k();
399  double lambda_r = (t - it.get_left()) / it.get_k();
400 
401  //Here we assume that the numbering of dofs goes from left to right!
402  result = *local_vectors[0];
403 
404  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
405  }
406  virtual void
407  InterpolateState(VECTOR& result,
408  const std::vector<VECTOR*> & local_vectors, double t,
409  const TimeIterator& it) const
410  {
411  assert(it.get_left() <= t);
412  assert(it.get_right() >= t);
413  if (local_vectors.size() != 2)
414  {
415  throw DOpEException(
416  "This function is currently not implemented for anything other than"
417  " linear interpolation of 2 DoFs.",
418  "MethodOfLine_SpaceTimeHandler::InterpolateState");
419  }
420  double lambda_l = (it.get_right() - t) / it.get_k();
421  double lambda_r = (t - it.get_left()) / it.get_k();
422 
423  //Here we assume that the numbering of dofs goes from left to right!
424  result = *local_vectors[0];
425 
426  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
427  }
428 
432  unsigned int
433  GetControlNDoFs(int /*time_point*/= -1) const
434  {
435  return GetControlDoFHandler().n_dofs();
436  }
440  unsigned int
441  GetStateNDoFs(int /*time_point*/= -1) const
442  {
443  return GetStateDoFHandler().n_dofs();
444  }
448  unsigned int
449  GetConstraintNDoFs(std::string name) const
450  {
451  return constraints_.n_dofs(name);
452  }
456  unsigned int
458  {
459  return constraints_.n_dofs("global");
460  //return constraints_.global_n_dofs();
461  }
465  unsigned int
467  {
468  //return constraints_.local_n_dofs();
469  return constraints_.n_dofs("local");
470  }
471 
475  const std::vector<Point<dealdim> >&
477  {
478  support_points_.resize(GetStateNDoFs());
480  GetStateDoFHandler(), support_points_);
481  return support_points_;
482  }
483 
484  /******************************************************/
489  void
490  ComputeControlSparsityPattern(SPARSITYPATTERN & sparsity) const;
491 
492  /******************************************************/
497  void
498  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
499  {
500  this->GetSparsityMaker()->ComputeSparsityPattern(
501  this->GetStateDoFHandler(), sparsity,
503  }
504 
505  /******************************************************/
506 
510  const FE<dealdim, dealdim>&
511  GetFESystem(std::string name) const
512  {
513  if (name == "state")
514  {
515  return *state_fe_;
516  }
517  else if (name == "control")
518  {
519  return *control_fe_;
520  }
521  else
522  {
523  throw DOpEException("Not implemented for name =" + name,
524  "MethodOfLines_SpaceTimeHandler::GetFESystem");
525  }
526 
527  }
528 
529  /******************************************************/
539  void
542  {
543  assert(ref_type == DOpEtypes::RefinementType::global);
544  RefineSpace(ref_type);
546  }
547 
548  /******************************************************/
559  void
562  {
563  assert(ref_type == DOpEtypes::RefinementType::global);
564  RefinementContainer ref_con_dummy;
565  RefineSpace(ref_con_dummy);
566  }
567 
568  /******************************************************/
578  void
579  RefineSpace(const RefinementContainer& ref_container)
580  {
581  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
582 
583  //make sure that we do not use any coarsening
584  assert(!ref_container.UsesCoarsening());
585 
586  if (control_mesh_transfer_ != NULL)
587  {
588  delete control_mesh_transfer_;
589  control_mesh_transfer_ = NULL;
590  }
591  if (state_mesh_transfer_ != NULL)
592  {
593  delete state_mesh_transfer_;
594  state_mesh_transfer_ = NULL;
595  }
596 #if dope_dimension == deal_II_dimension
597  control_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dopedim, VECTOR,
598  DH>(control_dof_handler_);
599 #endif
600  state_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dealdim, VECTOR,
601  DH>(state_dof_handler_);
602  if (DOpEtypes::RefinementType::global == ref_type)
603  {
604  triangulation_.set_all_refine_flags();
605  }
606  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
607  {
608  GridRefinement::refine_and_coarsen_fixed_number(triangulation_,
609  ref_container.GetLocalErrorIndicators(),
610  ref_container.GetTopFraction(),
611  ref_container.GetBottomFraction());
612  }
613  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
614  {
615  GridRefinement::refine_and_coarsen_fixed_fraction(triangulation_,
616  ref_container.GetLocalErrorIndicators(),
617  ref_container.GetTopFraction(),
618  ref_container.GetBottomFraction());
619  }
620  else if (DOpEtypes::RefinementType::optimized == ref_type)
621  {
622  GridRefinement::refine_and_coarsen_optimize(triangulation_,
623  ref_container.GetLocalErrorIndicators(),
624  ref_container.GetConvergenceOrder());
625  }
626  else
627  {
628  throw DOpEException("Not implemented for name =" + ref_type,
629  "MethodOfLines_SpaceTimeHandler::RefineSpace");
630  }
631  triangulation_.prepare_coarsening_and_refinement();
632 
633  //FIXME: works only if no coarsening happens, because we do not have the vectors to be interpolated availiable...
634  if (control_mesh_transfer_ != NULL)
635  control_mesh_transfer_->prepare_for_pure_refinement();
636  if (state_mesh_transfer_ != NULL)
637  state_mesh_transfer_->prepare_for_pure_refinement();
638 
639  triangulation_.execute_coarsening_and_refinement();
640  }
641 
642  /******************************************************/
643 
647  unsigned int
648  NewTimePointToOldTimePoint(unsigned int t) const
649  {
650  //TODO this has to be implemented when temporal refinement is possible!
651  //At present the temporal grid can't be refined
652  return t;
653  }
654 
655  /******************************************************/
656 
660  void
661  SpatialMeshTransferControl(const VECTOR& old_values,
662  VECTOR& new_values) const
663  {
664  if (control_mesh_transfer_ != NULL)
665  control_mesh_transfer_->refine_interpolate(old_values, new_values);
666  }
667  void
668  SpatialMeshTransferState(const VECTOR& old_values,
669  VECTOR& new_values) const
670  {
671  if (state_mesh_transfer_ != NULL)
672  state_mesh_transfer_->refine_interpolate(old_values, new_values);
673  }
674  /******************************************************/
683  void
686  {
687  user_defined_dof_constr_ = &constraints_maker;
688  user_defined_dof_constr_->RegisterMapping(this->GetMapping());
689  }
690  /******************************************************/
698  void
700  {
701  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_)
702  delete sparsitymaker_;
703  sparsitymaker_ = &sparsity_maker;
704  sparse_mkr_dynamic_ = false;
705  }
706 
707 // /******************************************************/
708 // /**
709 // * Through this function one can reinitialize the
710 // * triangulation for the state variable to be a copy of the
711 // * given argument.
712 // */
713 //
714 // void
715 // ResetTriangulation(const dealii::Triangulation<dealdim>& tria)
716 // {
717 // state_dof_handler_.clear();
718 // triangulation_.clear();
719 // triangulation_.copytriangulation_(tria);
720 // state_dof_handler_.initialize(triangulation_, *state_fe_);
721 // this->IncrementControlTicket();
722 // this->IncrementStateTicket();
723 // if (control_mesh_transfer_ != NULL)
724 // delete control_mesh_transfer_;
725 // control_mesh_transfer_ = NULL;
726 // if (state_mesh_transfer_ != NULL)
727 // delete state_mesh_transfer_;
728 // state_mesh_transfer_ = NULL;
729 //
730 // }
731 
732  private:
734  GetSparsityMaker() const
735  {
736  return sparsitymaker_;
737  }
738  const UserDefinedDoFConstraints<DH, dopedim, dealdim>*
739  GetUserDefinedDoFConstraints() const
740  {
741  return user_defined_dof_constr_;
742  }
743  SparsityMaker<DH, dealdim>* sparsitymaker_;
744  UserDefinedDoFConstraints<DH, dopedim, dealdim>* user_defined_dof_constr_;
745 
746  dealii::Triangulation<dealdim>& triangulation_;
747  DOpEWrapper::DoFHandler<dopedim, DH> control_dof_handler_;
748  DOpEWrapper::DoFHandler<dealdim, DH> state_dof_handler_;
749 
750  std::vector<unsigned int> control_dofs_per_block_;
751  std::vector<unsigned int> state_dofs_per_block_;
752 
753  dealii::ConstraintMatrix control_dof_constraints_;
754  dealii::ConstraintMatrix state_dof_constraints_;
755 
756  const dealii::SmartPointer<const FE<dealdim, dealdim> > control_fe_;
757  const dealii::SmartPointer<const FE<dealdim, dealdim> > state_fe_;
758 
759  const dealii::SmartPointer<const DOpEWrapper::Mapping<dealdim, DH> > mapping_;
760 
761  std::vector<Point<dealdim> > support_points_;
762 
763  Constraints constraints_;
766  bool sparse_mkr_dynamic_;
767  };
768 
769  /**************************explicit instantiation*************/
770 
774  template<>
775  void
776  DOpE::MethodOfLines_SpaceTimeHandler<dealii::FESystem,
777  dealii::DoFHandler, dealii::BlockSparsityPattern,
778  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
779  dealii::BlockSparsityPattern & sparsity) const;
780 
781  /******************************************************/
782 
783  template<>
784  void
785  MethodOfLines_SpaceTimeHandler<dealii::FESystem,
786  dealii::DoFHandler, dealii::SparsityPattern,
787  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
788  dealii::SparsityPattern & sparsity) const;
789 
790 
794  template<>
795  void
797  dealii::hp::FECollection,
798  dealii::hp::DoFHandler, dealii::BlockSparsityPattern,
799  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
800  dealii::BlockSparsityPattern & sparsity) const;
801 
802  /******************************************************/
803 
804  template<>
805  void
806  MethodOfLines_SpaceTimeHandler<dealii::hp::FECollection,
807  dealii::hp::DoFHandler, dealii::SparsityPattern,
808  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
809  dealii::SparsityPattern & sparsity) const;
810 
811 }
812 
813 #endif
const DOpEWrapper::DoFHandler< dopedim, DH > & GetControlDoFHandler() const
Definition: mol_spacetimehandler.h:293
void ReInitTime()
Definition: spacetimehandler_base.h:86
const DOpEWrapper::Mapping< dealdim, DH > & GetMapping() const
Definition: mol_spacetimehandler.h:311
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: spacetimehandler.h:327
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_spacetimehandler.h:498
unsigned int GetConstraintNDoFs(std::string name) const
Definition: mol_spacetimehandler.h:449
Definition: dopetypes.h:52
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_spacetimehandler.h:579
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_spacetimehandler.h:441
void SetSparsityMaker(SparsityMaker< DH, dealdim > &sparsity_maker)
Definition: mol_spacetimehandler.h:699
unsigned int GetStateDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_spacetimehandler.h:328
void SetActiveFEIndicesControl(DOpEWrapper::DoFHandler< dopedim, DH > &dof_handler)
Definition: spacetimehandler.h:348
void RefineSpaceTime(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_spacetimehandler.h:540
void ComputeControlSparsityPattern(SPARSITYPATTERN &sparsity) const
const std::vector< unsigned int > & GetConstraintDoFsPerBlock(std::string name) const
Definition: mol_spacetimehandler.h:360
unsigned int GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
Definition: mol_spacetimehandler.h:336
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:65
unsigned int GetControlDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_spacetimehandler.h:320
unsigned int GetControlNDoFs(int=-1) const
Definition: mol_spacetimehandler.h:433
const FE< dealdim, dealdim > & GetFESystem(std::string name) const
Definition: mol_spacetimehandler.h:511
ControlType
Definition: dopetypes.h:103
void SpatialMeshTransferControl(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_spacetimehandler.h:661
void SetUserDefinedDoFConstraints(UserDefinedDoFConstraints< DH, dopedim, dealdim > &constraints_maker)
Definition: mol_spacetimehandler.h:684
Definition: dopetypes.h:52
const std::vector< unsigned int > & GetControlDoFsPerBlock(int=-1) const
Definition: mol_spacetimehandler.h:344
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:98
const DOpEWrapper::DoFHandler< dealdim, DH > & GetStateDoFHandler() const
Definition: mol_spacetimehandler.h:302
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_spacetimehandler.h:648
Definition: timeiterator.h:63
Definition: userdefineddofconstraints.h:55
double get_k() const
Definition: timeiterator.h:195
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:95
RefinementType
Definition: dopetypes.h:50
virtual void InterpolateState(VECTOR &result, const std::vector< VECTOR * > &local_vectors, double t, const TimeIterator &it) const
Definition: mol_spacetimehandler.h:407
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, dealii::Triangulation< 1 > &times, const Constraints &c, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:163
Definition: refinementcontainer.h:42
const std::vector< unsigned int > & GetDoFsPerBlock(std::string name) const
Definition: constraints.h:159
virtual void InterpolateControl(VECTOR &result, const std::vector< VECTOR * > &local_vectors, double t, const TimeIterator &it) const
Definition: mol_spacetimehandler.h:385
Definition: active_fe_index_setter_interface.h:39
virtual const dealii::Vector< float > & GetLocalErrorIndicators() const
Definition: refinementcontainer.cc:46
virtual double GetTopFraction() const
Definition: refinementcontainer.cc:56
Definition: dopetypes.h:52
Definition: mol_spacetimehandler.h:51
Definition: constraints.h:33
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_spacetimehandler.h:352
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, const Constraints &c, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:130
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:85
double get_right() const
Definition: timeiterator.h:185
Definition: spacetimehandler.h:71
~MethodOfLines_SpaceTimeHandler()
Definition: mol_spacetimehandler.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 ReInit(std::vector< unsigned int > &control_dofs_per_block)
Definition: constraints.h:108
const dealii::ConstraintMatrix & GetControlDoFConstraints() const
Definition: mol_spacetimehandler.h:368
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
void ReInit(unsigned int control_n_blocks, const std::vector< unsigned int > &control_block_component, unsigned int state_n_blocks, const std::vector< unsigned int > &state_block_component)
Definition: mol_spacetimehandler.h:209
void RefineSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_spacetimehandler.h:560
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_spacetimehandler.h:668
bool UsesCoarsening() const
Definition: refinementcontainer.cc:103
const std::vector< Point< dealdim > > & GetMapDoFToSupportPoints()
Definition: mol_spacetimehandler.h:476
double get_left() const
Definition: timeiterator.h:175
void RefineTime(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: spacetimehandler_base.h:390
unsigned int GetNLocalConstraints() const
Definition: mol_spacetimehandler.h:466
void IncrementControlTicket()
Definition: spacetimehandler_base.h:446
Definition: solutiontransfer_wrapper.h:51
Definition: dopeexception.h:35
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_spacetimehandler.h:376
Definition: dopetypes.h:52
virtual double GetBottomFraction() const
Definition: refinementcontainer.cc:66
unsigned int n_dofs(std::string name) const
Definition: constraints.h:141
unsigned int GetNGlobalConstraints() const
Definition: mol_spacetimehandler.h:457
void IncrementStateTicket()
Definition: spacetimehandler_base.h:437