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 <deal.II/dofs/dof_handler.h>
37 #include <deal.II/dofs/dof_renumbering.h>
38 #include <deal.II/dofs/dof_tools.h>
39 #include <deal.II/hp/mapping_collection.h>
40 #include <deal.II/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:
66  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
67  const FE<dealdim, dealdim>& control_fe,
68  const FE<dealdim, dealdim>& state_fe,
70  bool flux_pattern = false,
73  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(type, index_setter),
74  triangulation_(triangulation),
75  control_dof_handler_(triangulation_),
76  state_dof_handler_(triangulation_),
77  control_fe_(&control_fe),
78  state_fe_(&state_fe),
79  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
80  constraints_(),
81  control_mesh_transfer_(NULL),
82  state_mesh_transfer_(NULL),
83  sparse_mkr_dynamic_(true)
84  {
85  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
86  user_defined_dof_constr_ = NULL;
87  }
88 
101  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
102  const FE<dealdim, dealdim>& control_fe,
103  const FE<dealdim, dealdim>& state_fe,
104  dealii::Triangulation<1> & times,
106  bool flux_pattern = false,
109  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(times, type, index_setter),
110  triangulation_(triangulation),
111  control_dof_handler_(triangulation_),
112  state_dof_handler_(triangulation_),
113  control_fe_(&control_fe),
114  state_fe_(&state_fe),
115  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
116  constraints_(),
117  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
118  {
119  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
120  user_defined_dof_constr_ = NULL;
121  }
122 
135  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
136  const FE<dealdim, dealdim>& control_fe,
137  const FE<dealdim, dealdim>& state_fe,
138  const Constraints& c,
140  bool flux_pattern = false,
143  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(type, index_setter),
144  triangulation_(triangulation),
145  control_dof_handler_(triangulation_),
146  state_dof_handler_(triangulation_),
147  control_fe_(&control_fe),
148  state_fe_(&state_fe),
149  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
150  constraints_(c),
151  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
152  {
153  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
154  user_defined_dof_constr_ = NULL;
155  }
156 
170  MethodOfLines_SpaceTimeHandler(dealii::Triangulation<dealdim>& triangulation,
171  const FE<dealdim, dealdim>& control_fe,
172  const FE<dealdim, dealdim>& state_fe,
173  dealii::Triangulation<1> & times,
174  const Constraints& c,
176  bool flux_pattern = false,
179  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim>(times, type, index_setter),
180  triangulation_(triangulation),
181  control_dof_handler_(triangulation_),
182  state_dof_handler_(triangulation_),
183  control_fe_(&control_fe),
184  state_fe_(&state_fe),
185  mapping_(&DOpEWrapper::StaticMappingQ1<dealdim, DH>::mapping_q1),
186  constraints_(c),
187  control_mesh_transfer_(NULL), state_mesh_transfer_(NULL), sparse_mkr_dynamic_(true)
188  {
189  sparsitymaker_ = new SparsityMaker<DH, dealdim>(flux_pattern);
190  user_defined_dof_constr_ = NULL;
191  }
192 
194  {
195  control_dof_handler_.clear();
196 
197  state_dof_handler_.clear();
198 
199  if (control_mesh_transfer_ != NULL)
200  {
201  delete control_mesh_transfer_;
202  }
203  if (state_mesh_transfer_ != NULL)
204  {
205  delete state_mesh_transfer_;
206  }
207  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_ == true)
208  {
209  delete sparsitymaker_;
210  }
211  }
212 
216  void
217  ReInit(unsigned int control_n_blocks,
218  const std::vector<unsigned int>& control_block_component,
219  unsigned int state_n_blocks,
220  const std::vector<unsigned int>& state_block_component)
221  {
222 #if dope_dimension > 0
223  SpaceTimeHandler<FE, DH, SPARSITYPATTERN,
224  VECTOR, dopedim, dealdim>::SetActiveFEIndicesControl(control_dof_handler_);
225 #endif
226  control_dof_handler_.distribute_dofs(*control_fe_);
227 
228 #if dope_dimension > 0
229  DoFRenumbering::component_wise (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_));
230  if(dopedim==dealdim)
231  {
232  control_dof_constraints_.clear ();
233  DoFTools::make_hanging_node_constraints (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_),
234  control_dof_constraints_);
235  if (GetUserDefinedDoFConstraints() != NULL)
236  GetUserDefinedDoFConstraints()->MakeControlDoFConstraints(control_dof_handler_,
237  control_dof_constraints_);
238  control_dof_constraints_.close ();
239  }
240  else
241  {
242  throw DOpEException("Not implemented for dopedim != dealdim","MethodOfLines_SpaceTimeHandler::ReInit");
243  }
244 #endif
245  control_dofs_per_block_.resize(control_n_blocks);
246 #if dope_dimension > 0
247  {
248  DoFTools::count_dofs_per_block (static_cast<DH<dopedim, dopedim>&>(control_dof_handler_),
249  control_dofs_per_block_,control_block_component);
250  }
251 #else
252  {
253  for (unsigned int i = 0; i < control_dofs_per_block_.size(); i++)
254  {
255  control_dofs_per_block_[i] = 0;
256  }
257  for (unsigned int i = 0; i < control_block_component.size(); i++)
258  {
259  control_dofs_per_block_[control_block_component[i]]++;
260  }
261  }
262 #endif
264  state_dof_handler_);
265  state_dof_handler_.distribute_dofs(GetFESystem("state"));
266  DoFRenumbering::component_wise(
267  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_));
268 
269  state_dof_constraints_.clear();
270  DoFTools::make_hanging_node_constraints(
271  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
272  state_dof_constraints_);
273  //TODO Dirichlet ueber Constraints
274  if (GetUserDefinedDoFConstraints() != NULL)
275  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
276  state_dof_handler_, state_dof_constraints_);
277  state_dof_constraints_.close();
278 
279  state_dofs_per_block_.resize(state_n_blocks);
280  DoFTools::count_dofs_per_block(
281  static_cast<DH<dealdim, dealdim>&>(state_dof_handler_),
282  state_dofs_per_block_, state_block_component);
283 
284  support_points_.clear();
285 
286  constraints_.ReInit(control_dofs_per_block_);
287  //constraints_.ReInit(control_dofs_per_block_, state_dofs_per_block_);
288 
289  //Initialize also the timediscretization.
290  this->ReInitTime();
291 
292  //There where changes invalidate tickets
293  this->IncrementControlTicket();
294  this->IncrementStateTicket();
295  }
296 
302  {
303  //There is only one mesh, hence always return this
304  return control_dof_handler_;
305  }
311  {
312  //There is only one mesh, hence always return this
313  return state_dof_handler_;
314  }
319  GetMapping() const
320  {
321  return *mapping_;
322  }
323 
327  unsigned int
328  GetControlDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
329  {
330  return control_dofs_per_block_[b];
331  }
335  unsigned int
336  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
337  {
338  return state_dofs_per_block_[b];
339  }
343  unsigned int
344  GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
345  {
346  return (constraints_.GetDoFsPerBlock(name))[b];
347  }
351  const std::vector<unsigned int>&
352  GetControlDoFsPerBlock(int /*time_point*/= -1) const
353  {
354  return control_dofs_per_block_;
355  }
359  const std::vector<unsigned int>&
360  GetStateDoFsPerBlock(int /*time_point*/= -1) const
361  {
362  return state_dofs_per_block_;
363  }
367  const std::vector<unsigned int>&
368  GetConstraintDoFsPerBlock(std::string name) const
369  {
370  return constraints_.GetDoFsPerBlock(name);
371  }
375  const dealii::ConstraintMatrix&
377  {
378  return control_dof_constraints_;
379  }
383  const dealii::ConstraintMatrix&
385  {
386  return state_dof_constraints_;
387  }
388 
392  virtual void
393  InterpolateControl(VECTOR& result,
394  const std::vector<VECTOR*> & local_vectors, double t,
395  const TimeIterator& it) const
396  {
397  assert(it.get_left() <= t);
398  assert(it.get_right() >= t);
399  if (local_vectors.size() != 2)
400  {
401  throw DOpEException(
402  "This function is currently not implemented for anything other than"
403  " linear interpolation of 2 DoFs.",
404  "MethodOfLine_SpaceTimeHandler::InterpolateControl");
405  }
406  double lambda_l = (it.get_right() - t) / it.get_k();
407  double lambda_r = (t - it.get_left()) / it.get_k();
408 
409  //Here we assume that the numbering of dofs goes from left to right!
410  result = *local_vectors[0];
411 
412  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
413  }
414  virtual void
415  InterpolateState(VECTOR& result,
416  const std::vector<VECTOR*> & local_vectors, double t,
417  const TimeIterator& it) const
418  {
419  assert(it.get_left() <= t);
420  assert(it.get_right() >= t);
421  if (local_vectors.size() != 2)
422  {
423  throw DOpEException(
424  "This function is currently not implemented for anything other than"
425  " linear interpolation of 2 DoFs.",
426  "MethodOfLine_SpaceTimeHandler::InterpolateState");
427  }
428  double lambda_l = (it.get_right() - t) / it.get_k();
429  double lambda_r = (t - it.get_left()) / it.get_k();
430 
431  //Here we assume that the numbering of dofs goes from left to right!
432  result = *local_vectors[0];
433 
434  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
435  }
436 
440  unsigned int
441  GetControlNDoFs(int /*time_point*/= -1) const
442  {
443  return GetControlDoFHandler().n_dofs();
444  }
448  unsigned int
449  GetStateNDoFs(int /*time_point*/= -1) const
450  {
451  return GetStateDoFHandler().n_dofs();
452  }
456  unsigned int
457  GetConstraintNDoFs(std::string name) const
458  {
459  return constraints_.n_dofs(name);
460  }
464  unsigned int
466  {
467  return constraints_.n_dofs("global");
468  //return constraints_.global_n_dofs();
469  }
473  unsigned int
475  {
476  //return constraints_.local_n_dofs();
477  return constraints_.n_dofs("local");
478  }
479 
483  const std::vector<Point<dealdim> >&
485  {
486  support_points_.resize(GetStateNDoFs());
488  GetStateDoFHandler(), support_points_);
489  return support_points_;
490  }
491 
492  /******************************************************/
497  void
498  ComputeControlSparsityPattern(SPARSITYPATTERN & sparsity) const;
499 
500  /******************************************************/
505  void
506  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
507  {
508  this->GetSparsityMaker()->ComputeSparsityPattern(
509  this->GetStateDoFHandler(), sparsity,
511  }
512 
513  /******************************************************/
514 
518  const FE<dealdim, dealdim>&
519  GetFESystem(std::string name) const
520  {
521  if (name == "state")
522  {
523  return *state_fe_;
524  }
525  else if (name == "control")
526  {
527  return *control_fe_;
528  }
529  else
530  {
531  throw DOpEException("Not implemented for name =" + name,
532  "MethodOfLines_SpaceTimeHandler::GetFESystem");
533  }
534 
535  }
536 
537  /******************************************************/
547  void
550  {
551  assert(ref_type == DOpEtypes::RefinementType::global);
552  RefineSpace(ref_type);
554  }
555 
556  /******************************************************/
567  void
570  {
571  //assert(ref_type == DOpEtypes::RefinementType::global);
572  RefinementContainer ref_con_dummy;
573  RefineSpace(ref_con_dummy);
574  }
575 
576  /******************************************************/
586  void
587  RefineSpace(const RefinementContainer& ref_container)
588  {
589  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
590 
591  //make sure that we do not use any coarsening
592  assert(!ref_container.UsesCoarsening());
593 
594  if (control_mesh_transfer_ != NULL)
595  {
596  delete control_mesh_transfer_;
597  control_mesh_transfer_ = NULL;
598  }
599  if (state_mesh_transfer_ != NULL)
600  {
601  delete state_mesh_transfer_;
602  state_mesh_transfer_ = NULL;
603  }
604 #if dope_dimension == deal_II_dimension
605  control_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dopedim, VECTOR,
606  DH>(control_dof_handler_);
607 #endif
608  state_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dealdim, VECTOR,
609  DH>(state_dof_handler_);
610  if (DOpEtypes::RefinementType::global == ref_type)
611  {
612  triangulation_.set_all_refine_flags();
613  }
614  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
615  {
616  GridRefinement::refine_and_coarsen_fixed_number(triangulation_,
617  ref_container.GetLocalErrorIndicators(),
618  ref_container.GetTopFraction(),
619  ref_container.GetBottomFraction());
620  }
621  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
622  {
623  GridRefinement::refine_and_coarsen_fixed_fraction(triangulation_,
624  ref_container.GetLocalErrorIndicators(),
625  ref_container.GetTopFraction(),
626  ref_container.GetBottomFraction());
627  }
628  else if (DOpEtypes::RefinementType::optimized == ref_type)
629  {
630  GridRefinement::refine_and_coarsen_optimize(triangulation_,
631  ref_container.GetLocalErrorIndicators(),
632  ref_container.GetConvergenceOrder());
633  }
634  else
635  {
636  throw DOpEException("Not implemented for name =" + DOpEtypesToString(ref_type),
637  "MethodOfLines_SpaceTimeHandler::RefineSpace");
638  }
639  triangulation_.prepare_coarsening_and_refinement();
640 
641  //FIXME: works only if no coarsening happens, because we do not have the vectors to be interpolated availiable...
642  if (control_mesh_transfer_ != NULL)
643  control_mesh_transfer_->prepare_for_pure_refinement();
644  if (state_mesh_transfer_ != NULL)
645  state_mesh_transfer_->prepare_for_pure_refinement();
646 
647  triangulation_.execute_coarsening_and_refinement();
648  }
649 
650  /******************************************************/
651 
655  unsigned int
656  NewTimePointToOldTimePoint(unsigned int t) const
657  {
658  //TODO this has to be implemented when temporal refinement is possible!
659  //At present the temporal grid can't be refined
660  return t;
661  }
662 
663  /******************************************************/
664 
668  void
669  SpatialMeshTransferControl(const VECTOR& old_values,
670  VECTOR& new_values) const
671  {
672  if (control_mesh_transfer_ != NULL)
673  control_mesh_transfer_->refine_interpolate(old_values, new_values);
674  }
675  void
676  SpatialMeshTransferState(const VECTOR& old_values,
677  VECTOR& new_values) const
678  {
679  if (state_mesh_transfer_ != NULL)
680  state_mesh_transfer_->refine_interpolate(old_values, new_values);
681  }
682  /******************************************************/
691  void
694  {
695  user_defined_dof_constr_ = &constraints_maker;
696  user_defined_dof_constr_->RegisterMapping(this->GetMapping());
697  }
698  /******************************************************/
706  void
708  {
709  assert(sparse_mkr_dynamic_==true); //If not true, we already set the sparsity maker
710  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_)
711  delete sparsitymaker_;
712  sparsitymaker_ = &sparsity_maker;
713  sparse_mkr_dynamic_ = false;
714  }
715 
716  /******************************************************/
722  void
723  ResetTriangulation(const dealii::Triangulation<dealdim>& tria);
724 
725  private:
727  GetSparsityMaker() const
728  {
729  return sparsitymaker_;
730  }
731  const UserDefinedDoFConstraints<DH, dopedim, dealdim>*
732  GetUserDefinedDoFConstraints() const
733  {
734  return user_defined_dof_constr_;
735  }
736  SparsityMaker<DH, dealdim>* sparsitymaker_;
737  UserDefinedDoFConstraints<DH, dopedim, dealdim>* user_defined_dof_constr_;
738 
739  dealii::Triangulation<dealdim>& triangulation_;
740  DOpEWrapper::DoFHandler<dopedim, DH> control_dof_handler_;
741  DOpEWrapper::DoFHandler<dealdim, DH> state_dof_handler_;
742 
743  std::vector<unsigned int> control_dofs_per_block_;
744  std::vector<unsigned int> state_dofs_per_block_;
745 
746  dealii::ConstraintMatrix control_dof_constraints_;
747  dealii::ConstraintMatrix state_dof_constraints_;
748 
749  const dealii::SmartPointer<const FE<dealdim, dealdim> > control_fe_;
750  const dealii::SmartPointer<const FE<dealdim, dealdim> > state_fe_;
751 
752  const dealii::SmartPointer<const DOpEWrapper::Mapping<dealdim, DH> > mapping_;
753 
754  std::vector<Point<dealdim> > support_points_;
755 
756  Constraints constraints_;
759  bool sparse_mkr_dynamic_;
760  };
761 
762  /**************************explicit instantiation*************/
763 
767  template<>
768  void
769  DOpE::MethodOfLines_SpaceTimeHandler<dealii::FESystem,
770  dealii::DoFHandler, dealii::BlockSparsityPattern,
771  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
772  dealii::BlockSparsityPattern & sparsity) const;
773  template<>
774  void
775  DOpE::MethodOfLines_SpaceTimeHandler<dealii::FESystem,
776  dealii::DoFHandler, dealii::BlockSparsityPattern,
777  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ResetTriangulation(
778  const dealii::Triangulation<deal_II_dimension>& tria);
779 
780  /******************************************************/
781 
782  template<>
783  void
784  MethodOfLines_SpaceTimeHandler<dealii::FESystem,
785  dealii::DoFHandler, dealii::SparsityPattern,
786  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
787  dealii::SparsityPattern & sparsity) const;
788  template<>
789  void
790  MethodOfLines_SpaceTimeHandler<dealii::FESystem,
791  dealii::DoFHandler, dealii::SparsityPattern,
792  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ResetTriangulation(
793  const dealii::Triangulation<deal_II_dimension>& tria);
794 
795 
799  template<>
800  void
802  dealii::hp::FECollection,
803  dealii::hp::DoFHandler, dealii::BlockSparsityPattern,
804  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
805  dealii::BlockSparsityPattern & sparsity) const;
806  template<>
807  void
809  dealii::hp::FECollection,
810  dealii::hp::DoFHandler, dealii::BlockSparsityPattern,
811  dealii::BlockVector<double>, dope_dimension, deal_II_dimension>::ResetTriangulation(
812  const dealii::Triangulation<deal_II_dimension>& tria);
813 
814  /******************************************************/
815 
816  template<>
817  void
818  MethodOfLines_SpaceTimeHandler<dealii::hp::FECollection,
819  dealii::hp::DoFHandler, dealii::SparsityPattern,
820  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ComputeControlSparsityPattern(
821  dealii::SparsityPattern & sparsity) const;
822  template<>
823  void
824  MethodOfLines_SpaceTimeHandler<dealii::hp::FECollection,
825  dealii::hp::DoFHandler, dealii::SparsityPattern,
826  dealii::Vector<double>, dope_dimension, deal_II_dimension>::ResetTriangulation(
827  const dealii::Triangulation<deal_II_dimension>& tria);
828 
829 }
830 
831 #endif
const DOpEWrapper::DoFHandler< dopedim, DH > & GetControlDoFHandler() const
Definition: mol_spacetimehandler.h:301
void ReInitTime()
Definition: spacetimehandler_base.h:86
const DOpEWrapper::Mapping< dealdim, DH > & GetMapping() const
Definition: mol_spacetimehandler.h:319
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: spacetimehandler.h:327
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_spacetimehandler.h:506
unsigned int GetConstraintNDoFs(std::string name) const
Definition: mol_spacetimehandler.h:457
Definition: dopetypes.h:52
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, const Constraints &c, DOpEtypes::ControlType type, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:135
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_spacetimehandler.h:587
void RefineTime(DOpEtypes::RefinementType=DOpEtypes::RefinementType::global)
Definition: spacetimehandler_base.h:390
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_spacetimehandler.h:449
void SetSparsityMaker(SparsityMaker< DH, dealdim > &sparsity_maker)
Definition: mol_spacetimehandler.h:707
unsigned int GetStateDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_spacetimehandler.h:336
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:548
void ComputeControlSparsityPattern(SPARSITYPATTERN &sparsity) const
const std::vector< unsigned int > & GetConstraintDoFsPerBlock(std::string name) const
Definition: mol_spacetimehandler.h:368
unsigned int GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
Definition: mol_spacetimehandler.h:344
unsigned int GetControlDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_spacetimehandler.h:328
unsigned int GetControlNDoFs(int=-1) const
Definition: mol_spacetimehandler.h:441
const FE< dealdim, dealdim > & GetFESystem(std::string name) const
Definition: mol_spacetimehandler.h:519
ControlType
Definition: dopetypes.h:103
void SpatialMeshTransferControl(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_spacetimehandler.h:669
void SetUserDefinedDoFConstraints(UserDefinedDoFConstraints< DH, dopedim, dealdim > &constraints_maker)
Definition: mol_spacetimehandler.h:692
Definition: dopetypes.h:52
const std::vector< unsigned int > & GetControlDoFsPerBlock(int=-1) const
Definition: mol_spacetimehandler.h:352
const DOpEWrapper::DoFHandler< dealdim, DH > & GetStateDoFHandler() const
Definition: mol_spacetimehandler.h:310
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_spacetimehandler.h:656
Definition: timeiterator.h:62
Definition: userdefineddofconstraints.h:55
double get_k() const
Definition: timeiterator.h:194
void ResetTriangulation(const dealii::Triangulation< dealdim > &tria)
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:86
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:415
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, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:101
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:393
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
void RefineSpace(DOpEtypes::RefinementType=DOpEtypes::RefinementType::global)
Definition: mol_spacetimehandler.h:568
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_spacetimehandler.h:360
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:76
double get_right() const
Definition: timeiterator.h:184
Definition: spacetimehandler.h:71
~MethodOfLines_SpaceTimeHandler()
Definition: mol_spacetimehandler.h:193
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:376
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
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:217
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, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:170
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_spacetimehandler.h:676
bool UsesCoarsening() const
Definition: refinementcontainer.cc:94
const std::vector< Point< dealdim > > & GetMapDoFToSupportPoints()
Definition: mol_spacetimehandler.h:484
double get_left() const
Definition: timeiterator.h:174
unsigned int GetNLocalConstraints() const
Definition: mol_spacetimehandler.h:474
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:384
Definition: dopetypes.h:52
virtual double GetBottomFraction() const
Definition: refinementcontainer.cc:66
std::string DOpEtypesToString(const C &)
Definition: dopetypes.h:134
MethodOfLines_SpaceTimeHandler(dealii::Triangulation< dealdim > &triangulation, const FE< dealdim, dealdim > &control_fe, const FE< dealdim, dealdim > &state_fe, DOpEtypes::ControlType type, bool flux_pattern=false, const ActiveFEIndexSetterInterface< dopedim, dealdim > &index_setter=ActiveFEIndexSetterInterface< dopedim, dealdim >())
Definition: mol_spacetimehandler.h:66
unsigned int n_dofs(std::string name) const
Definition: constraints.h:141
unsigned int GetNGlobalConstraints() const
Definition: mol_spacetimehandler.h:465
void IncrementStateTicket()
Definition: spacetimehandler_base.h:437