DOpE
mol_multimesh_spacetimehandler.h
Go to the documentation of this file.
1 
24 #ifndef _MOL_MULTIMESH_SPACE_TIME_HANDLER_H_
25 #define _MOL_MULTIMESH_SPACE_TIME_HANDLER_H_
26 
27 #include "mol_spacetimehandler.h"
28 #include "constraints.h"
29 #include "sparsitymaker.h"
31 #include "sth_internals.h"
32 #include "refinementcontainer.h"
34 
35 #include <dofs/dof_handler.h>
36 #include <dofs/dof_renumbering.h>
37 #include <dofs/dof_tools.h>
38 #include <hp/mapping_collection.h>
39 #include <lac/constraint_matrix.h>
40 #include <deal.II/grid/grid_refinement.h>
41 
42 namespace DOpE
43 {
55  template<template<int, int> class FE, template<int, int> class DH,
56  typename SPARSITYPATTERN, typename VECTOR, int dim>
58  DH, SPARSITYPATTERN, VECTOR, dim, dim>
59  {
60  public:
62  dealii::Triangulation<dim>& triangulation,
63  const FE<dim, dim>& control_fe, const FE<dim, dim>& state_fe,
65  const ActiveFEIndexSetterInterface<dim, dim>& index_setter =
67  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dim, dim>(type,
68  index_setter), _state_triangulation(triangulation), _control_dof_handler(
69  _control_triangulation), _state_dof_handler(
70  _state_triangulation), _control_fe(&control_fe), _state_fe(
71  &state_fe), _mapping(
72  DOpEWrapper::StaticMappingQ1<dim, DH>::mapping_q1), _constraints(), _control_mesh_transfer(
73  NULL), _state_mesh_transfer(NULL), _sparse_mkr_dynamic(true)
74  {
75  _control_triangulation.copy_triangulation(_state_triangulation);
76  _sparsitymaker = new SparsityMaker<DH, dim>;
77  _user_defined_dof_constr = NULL;
78  }
79 
81  dealii::Triangulation<dim>& triangulation,
82  const FE<dim, dim>& control_fe, const FE<dim, dim>& state_fe,
83  const dealii::Triangulation<1> & times, DOpEtypes::ControlType type,
84  const ActiveFEIndexSetterInterface<dim, dim>& index_setter =
86  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dim, dim>(times,
87  type, index_setter), _state_triangulation(triangulation), _control_dof_handler(
88  _control_triangulation), _state_dof_handler(
89  _state_triangulation), _control_fe(&control_fe), _state_fe(
90  &state_fe), _mapping(
91  DOpEWrapper::StaticMappingQ1<dim, DH>::mapping_q1), _constraints(), _control_mesh_transfer(
92  NULL), _state_mesh_transfer(NULL), _sparse_mkr_dynamic(true)
93  {
94  _control_triangulation.copy_triangulation(_state_triangulation);
95  _sparsitymaker = new SparsityMaker<DH, dim>;
96  _user_defined_dof_constr = NULL;
97  }
98 
100  dealii::Triangulation<dim>& triangulation,
101  const FE<dim, dim>& control_fe, const FE<dim, dim>& state_fe,
102  const Constraints& c, DOpEtypes::ControlType type,
103  const ActiveFEIndexSetterInterface<dim, dim>& index_setter =
105  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dim, dim>(type,
106  index_setter), _state_triangulation(triangulation), _control_dof_handler(
107  _control_triangulation), _state_dof_handler(
108  _state_triangulation), _control_fe(&control_fe), _state_fe(
109  &state_fe), _mapping(
110  DOpEWrapper::StaticMappingQ1<dim, DH>::mapping_q1), _constraints(
111  c), _control_mesh_transfer(NULL), _state_mesh_transfer(NULL), _sparse_mkr_dynamic(
112  true)
113  {
114  _control_triangulation.copy_triangulation(_state_triangulation);
115  _sparsitymaker = new SparsityMaker<DH, dim>;
116  _user_defined_dof_constr = NULL;
117  }
118 
120  dealii::Triangulation<dim>& triangulation,
121  const FE<dim, dim>& control_fe, const FE<dim, dim>& state_fe,
122  const dealii::Triangulation<1> & times, const Constraints& c,
124  const ActiveFEIndexSetterInterface<dim, dim>& index_setter =
126  SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR, dim, dim>(times,
127  type, index_setter), _state_triangulation(triangulation), _control_dof_handler(
128  _control_triangulation), _state_dof_handler(
129  _state_triangulation), _control_fe(&control_fe), _state_fe(
130  &state_fe), _mapping(
131  DOpEWrapper::StaticMappingQ1<dim, DH>::mapping_q1), _constraints(
132  c), _control_mesh_transfer(NULL), _state_mesh_transfer(NULL), _sparse_mkr_dynamic(
133  true)
134  {
135  _control_triangulation.copy_triangulation(_state_triangulation);
136  _sparsitymaker = new SparsityMaker<DH, dim>;
137  _user_defined_dof_constr = NULL;
138  }
139 
141  {
142  _control_dof_handler.clear();
143 
144  _state_dof_handler.clear();
145 
146  if (_control_mesh_transfer != NULL)
147  {
148  delete _control_mesh_transfer;
149  }
150  if (_state_mesh_transfer != NULL)
151  {
152  delete _state_mesh_transfer;
153  }
154  if (_sparsitymaker != NULL && _sparse_mkr_dynamic == true)
155  {
156  delete _sparsitymaker;
157  }
158  }
159 
163  void
164  ReInit(unsigned int control_n_blocks,
165  const std::vector<unsigned int>& control_block_component,
166  unsigned int state_n_blocks,
167  const std::vector<unsigned int>& state_block_component)
168  {
170  _control_dof_handler);
171  _control_dof_handler.distribute_dofs(*_control_fe);
172  DoFRenumbering::component_wise(static_cast<DH<dim, dim>&>(_control_dof_handler));
173 
174  _control_dof_constraints.clear();
175  DoFTools::make_hanging_node_constraints(
176  static_cast<DH<dim, dim>&>(_control_dof_handler), _control_dof_constraints);
177  if (GetUserDefinedDoFConstraints() != NULL)
178  GetUserDefinedDoFConstraints()->MakeControlDoFConstraints(
179  _control_dof_handler, _control_dof_constraints);
180  _control_dof_constraints.close();
181 
182  _control_dofs_per_block.resize(control_n_blocks);
183  {
184  DoFTools::count_dofs_per_block(
185  static_cast<DH<dim, dim>&>(_control_dof_handler), _control_dofs_per_block,
186  control_block_component);
187  }
188 
190  _state_dof_handler);
191  _state_dof_handler.distribute_dofs(GetFESystem("state"));
192  DoFRenumbering::component_wise(static_cast<DH<dim, dim>&>(_state_dof_handler));
193 
194  _state_dof_constraints.clear();
195  DoFTools::make_hanging_node_constraints(
196  static_cast<DH<dim, dim>&>(_state_dof_handler), _state_dof_constraints);
197  //TODO Dirichlet Daten hierueber.
198  if (GetUserDefinedDoFConstraints() != NULL)
199  GetUserDefinedDoFConstraints()->MakeStateDoFConstraints(
200  _state_dof_handler, _state_dof_constraints);
201  _state_dof_constraints.close();
202 
203  _state_dofs_per_block.resize(state_n_blocks);
204  DoFTools::count_dofs_per_block(static_cast<DH<dim, dim>&>(_state_dof_handler),
205  _state_dofs_per_block, state_block_component);
206 
207  _support_points.clear();
208 
209  _constraints.ReInit(_control_dofs_per_block);
210  //_constraints.ReInit(_control_dofs_per_block, _state_dofs_per_block);
211 
212  //Initialize also the timediscretization.
213  this->ReInitTime();
214 
215  //There where changes invalidate tickets
216  this->IncrementControlTicket();
217  this->IncrementStateTicket();
218  }
219 
225  {
226  //There is only one mesh, hence always return this
227  return _control_dof_handler;
228  }
234  {
235  //There is only one mesh, hence always return this
236  return _state_dof_handler;
237  }
238 
243  GetMapping() const
244  {
245  return _mapping;
246  }
247 
251  unsigned int
252  GetControlDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
253  {
254  return _control_dofs_per_block[b];
255  }
259  unsigned int
260  GetStateDoFsPerBlock(unsigned int b, int /*time_point*/= -1) const
261  {
262  return _state_dofs_per_block[b];
263  }
267  unsigned int
268  GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
269  {
270  return (_constraints.GetDoFsPerBlock(name))[b];
271  }
275  const std::vector<unsigned int>&
277  {
278  return _control_dofs_per_block;
279  }
283  const std::vector<unsigned int>&
284  GetStateDoFsPerBlock(int /*time_point*/= -1) const
285  {
286  return _state_dofs_per_block;
287  }
291  const std::vector<unsigned int>&
292  GetConstraintDoFsPerBlock(std::string name) const
293  {
294  return _constraints.GetDoFsPerBlock(name);
295  }
299  const dealii::ConstraintMatrix&
301  {
302  return _control_dof_constraints;
303  }
307  const dealii::ConstraintMatrix&
309  {
310  return _state_dof_constraints;
311  }
312 
317  virtual void
318  InterpolateState(VECTOR& result,
319  const std::vector<VECTOR*> & local_vectors, double t,
320  const TimeIterator& it) const
321  {
322  assert(it.get_left() <= t);
323  assert(it.get_right() >= t);
324  if (local_vectors.size() != 2)
325  throw DOpEException(
326  "This function is currently not implemented for anything other than"
327  " linear interpolation of 2 DoFs.",
328  "MethodOfLine_SpaceTimeHandler::InterpolateState");
329 
330  double lambda_l = (it.get_right() - t) / it.get_k();
331  double lambda_r = (t - it.get_left()) / it.get_k();
332 
333  //Here we assume that the numbering of dofs goes from left to right!
334  result = *local_vectors[0];
335 
336  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
337  }
338 
342  unsigned int
344  {
345  return GetControlDoFHandler().n_dofs();
346  }
350  unsigned int
351  GetStateNDoFs(int /*time_point*/= -1) const
352  {
353  return GetStateDoFHandler().n_dofs();
354  }
358  unsigned int
359  GetConstraintNDoFs(std::string name) const
360  {
361  return _constraints.n_dofs(name);
362  }
366  unsigned int
368  {
369  return _constraints.n_dofs("global");
370  //return _constraints.global_n_dofs();
371  }
375  unsigned int
377  {
378  return _constraints.n_dofs("local");
379  //return _constraints.local_n_dofs();
380  }
381 
385  const std::vector<Point<dim> >&
387  {
388  _support_points.resize(GetStateNDoFs());
390  GetStateDoFHandler(), _support_points);
391  return _support_points;
392  }
393 
394  /******************************************************/
395  void
396  ComputeControlSparsityPattern(SPARSITYPATTERN & sparsity) const;
397 
398  /******************************************************/
399  void
400  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
401  {
402  this->GetSparsityMaker()->ComputeSparsityPattern(
403  this->GetStateDoFHandler(), sparsity,
405  }
406 
407  /******************************************************/
408 
412  const FE<dim, dim>&
413  GetFESystem(std::string name) const
414  {
415  if (name == "state")
416  {
417  return *_state_fe;
418  }
419  else if (name == "control")
420  {
421  return *_control_fe;
422  }
423  else
424  {
425  throw DOpEException("Not implemented for name =" + name,
426  "MethodOfLines_MultiMesh_SpaceTimeHandler::GetFESystem");
427  }
428 
429  }
430 
441  void
444  {
445  assert(ref_type == DOpEtypes::RefinementType::global);
446  RefinementContainer ref_con_dummy;
447  RefineStateSpace(ref_con_dummy);
448  RefineControlSpace(ref_con_dummy);
449  }
450 
461  void
464  {
465  assert(ref_type == DOpEtypes::RefinementType::global);
466  RefinementContainer ref_con_dummy;
467  RefineStateSpace(ref_con_dummy);
468  }
469 
480  void
483  {
484  assert(ref_type == DOpEtypes::RefinementType::global);
485  RefinementContainer ref_con_dummy;
486  RefineControlSpace(ref_con_dummy);
487  }
488 
498  template<typename NUMBER>
499  void
500  RefineSpace(const RefinementContainer&ref_container)
501  {
502  RefineStateSpace(ref_container);
503  RefineControlSpace(ref_container);
504  }
505 
514  void
515  RefineStateSpace(const RefinementContainer& ref_container)
516  {
517  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
518 
519  //make sure that we do not use any coarsening
520  assert(!ref_container.UsesCoarsening());
521 
522  if (_state_mesh_transfer != NULL)
523  {
524  delete _state_mesh_transfer;
525  _state_mesh_transfer = NULL;
526  }
527  _state_mesh_transfer = new DOpEWrapper::SolutionTransfer<dim, VECTOR,
528  DH>(_state_dof_handler);
529 
530  if (DOpEtypes::RefinementType::global == ref_type)
531  {
532  _state_triangulation.set_all_refine_flags();
533  }
534  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
535  {
536 
537  GridRefinement::refine_and_coarsen_fixed_number(
538  _state_triangulation, ref_container.GetLocalErrorIndicators(),
539  ref_container.GetTopFraction(),
540  ref_container.GetBottomFraction());
541  }
542  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
543  {
544 
545  GridRefinement::refine_and_coarsen_fixed_fraction(
546  _state_triangulation, ref_container.GetLocalErrorIndicators(),
547  ref_container.GetTopFraction(),
548  ref_container.GetBottomFraction());
549  }
550  else if (DOpEtypes::RefinementType::optimized == ref_type)
551  {
552 
553  GridRefinement::refine_and_coarsen_optimize(_state_triangulation,
554  ref_container.GetLocalErrorIndicators(),
555  ref_container.GetConvergenceOrder());
556  }
557  else if (DOpEtypes::RefinementType::finest_of_both == ref_type)
558  {
559  this->FlagIfLeftIsNotFinest(_state_triangulation,
560  _control_triangulation);
561  }
562  else
563  {
564  throw DOpEException("Not implemented for name =" + ref_type,
565  "MethodOfLines_MultiMesh_SpaceTimeHandler::RefineSpace");
566  }
567  _state_triangulation.prepare_coarsening_and_refinement();
568  if (_state_mesh_transfer != NULL)
569  _state_mesh_transfer->prepare_for_pure_refinement();
570 
571  _state_triangulation.execute_coarsening_and_refinement();
572  }
573 
582  void
584  {
585  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
586 
587  //make sure that we do not use any coarsening
588  assert(!ref_container.UsesCoarsening());
589 
590  if (_control_mesh_transfer != NULL)
591  {
592  delete _control_mesh_transfer;
593  _control_mesh_transfer = NULL;
594  }
595 #if dope_dimension == deal_II_dimension
596  _control_mesh_transfer = new DOpEWrapper::SolutionTransfer<dim, VECTOR,
597  DH>(_control_dof_handler);
598 #endif
599  if (DOpEtypes::RefinementType::global == ref_type)
600  {
601  _control_triangulation.set_all_refine_flags();
602  }
603 
604  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
605  {
606 
607  GridRefinement::refine_and_coarsen_fixed_number(
608  _control_triangulation, ref_container.GetLocalErrorIndicators(),
609  ref_container.GetTopFraction(),
610  ref_container.GetBottomFraction());
611  }
612  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
613  {
614 
615  GridRefinement::refine_and_coarsen_fixed_fraction(
616  _control_triangulation, ref_container.GetLocalErrorIndicators(),
617  ref_container.GetTopFraction(),
618  ref_container.GetBottomFraction());
619  }
620  else if (DOpEtypes::RefinementType::optimized == ref_type)
621  {
622 
623  GridRefinement::refine_and_coarsen_optimize(_control_triangulation,
624  ref_container.GetLocalErrorIndicators(),
625  ref_container.GetConvergenceOrder());
626  }
627  else if ("finest-of-both")
628  {
629  this->FlagIfLeftIsNotFinest(_control_triangulation,
630  _state_triangulation);
631  }
632  else
633  {
634  throw DOpEException("Not implemented for name =" + ref_type,
635  "MethodOfLines_MultiMesh_SpaceTimeHandler::RefineSpace");
636  }
637  _control_triangulation.prepare_coarsening_and_refinement();
638 
639  //FIXME: works only if no coarsening happens, because we do not have the vectors to be interpolated availiable...
640  if (_control_mesh_transfer != NULL)
641  _control_mesh_transfer->prepare_for_pure_refinement();
642 
643  _control_triangulation.execute_coarsening_and_refinement();
644  }
645  /******************************************************/
646 
650  unsigned int
651  NewTimePointToOldTimePoint(unsigned int t) const
652  {
653  //TODO this has to be implemented when temporal refinement is possible!
654  //At present the temporal grid can't be refined
655  return t;
656  }
657 
658  /******************************************************/
659 
663  void
664  SpatialMeshTransferControl(const VECTOR& old_values,
665  VECTOR& new_values) const
666  {
667  if (_control_mesh_transfer != NULL)
668  _control_mesh_transfer->refine_interpolate(old_values, new_values);
669  }
670  /******************************************************/
671 
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
693  UserDefinedDoFConstraints<DH, dim>& constraints_maker)
694  {
695  _user_defined_dof_constr = &constraints_maker;
696  _user_defined_dof_constr->RegisterMapping(this->GetMapping());
697  }
698  /******************************************************/
706  void
708  {
709  if (_sparsitymaker != NULL && _sparse_mkr_dynamic)
710  delete _sparsitymaker;
711  _sparsitymaker = &sparsity_maker;
712  _sparse_mkr_dynamic = false;
713  }
714  /******************************************************/
721  void
722  ResetStateTriangulation(const dealii::Triangulation<dim>& tria)
723  {
724  _state_dof_handler.clear();
725  _state_triangulation.clear();
726  _state_triangulation.copy_triangulation(tria);
727  _state_dof_handler.initialize(_state_triangulation, *_state_fe);
728  this->IncrementStateTicket();
729  if (_state_mesh_transfer != NULL)
730  delete _state_mesh_transfer;
731  _state_mesh_transfer = NULL;
732  }
733 
734  /******************************************************/
741  void
742  ResetControlTriangulation(const dealii::Triangulation<dim>& tria)
743  {
744  _control_dof_handler.clear();
745  _control_triangulation.clear();
746  _control_triangulation.copy_triangulation(tria);
747  _control_dof_handler.initialize(_control_triangulation, *_control_fe);
748  this->IncrementControlTicket();
749  if (_control_mesh_transfer != NULL)
750  delete _control_mesh_transfer;
751  _control_mesh_transfer = NULL;
752  }
753 
754  private:
756  GetSparsityMaker() const
757  {
758  return _sparsitymaker;
759  }
760  const UserDefinedDoFConstraints<DH, dim>*
761  GetUserDefinedDoFConstraints() const
762  {
763  return _user_defined_dof_constr;
764  }
765 
766  void
767  FlagIfLeftIsNotFinest(dealii::Triangulation<dim>& left,
768  const dealii::Triangulation<dim>& right);
769 
770  SparsityMaker<DH, dim>* _sparsitymaker;
771  UserDefinedDoFConstraints<DH, dim>* _user_defined_dof_constr;
772 
773  dealii::Triangulation<dim> _control_triangulation;
774  dealii::Triangulation<dim>& _state_triangulation;
775  DOpEWrapper::DoFHandler<dim, DH> _control_dof_handler;
776  DOpEWrapper::DoFHandler<dim, DH> _state_dof_handler;
777 
778  std::vector<unsigned int> _control_dofs_per_block;
779  std::vector<unsigned int> _state_dofs_per_block;
780 
781  dealii::ConstraintMatrix _control_dof_constraints;
782  dealii::ConstraintMatrix _state_dof_constraints;
783 
784  const dealii::SmartPointer<const FE<dim, dim>> _control_fe;
785  const dealii::SmartPointer<const FE<dim, dim>> _state_fe;
786 
788 
789  std::vector<Point<dim> > _support_points;
790 
791  Constraints _constraints;
792  DOpEWrapper::SolutionTransfer<dim, VECTOR,DH>* _control_mesh_transfer;
793  DOpEWrapper::SolutionTransfer<dim, VECTOR,DH>* _state_mesh_transfer;
794  bool _sparse_mkr_dynamic;
795  };
796 
797  /**************************explicit instantiation*************/
798 #if dope_dimension == deal_II_dimension
799 
802  template<>
803  void
805  dealii::DoFHandler, dealii::BlockSparsityPattern,
806  dealii::BlockVector<double>, dope_dimension>::ComputeControlSparsityPattern(
807  dealii::BlockSparsityPattern & sparsity) const
808  {
809  const std::vector<unsigned int>& blocks = this->GetControlDoFsPerBlock();
810  dealii::BlockCompressedSimpleSparsityPattern csp(blocks.size(),
811  blocks.size());
812  for (unsigned int i = 0; i < blocks.size(); i++)
813  {
814  for (unsigned int j = 0; j < blocks.size(); j++)
815  {
816  csp.block(i, j).reinit(this->GetControlDoFsPerBlock(i),
817  this->GetControlDoFsPerBlock(j));
818  }
819  }
820  csp.collect_sizes();
821 
822  dealii::DoFTools::make_sparsity_pattern(
823  static_cast<const dealii::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
824  csp);
825  this->GetControlDoFConstraints().condense(csp);
826  sparsity.copy_from(csp);
827  }
828 
829  /******************************************************/
830 
831  template<>
832  void
834  dealii::DoFHandler, dealii::SparsityPattern, dealii::Vector<double>,
835  dope_dimension>::ComputeControlSparsityPattern(
836  dealii::SparsityPattern & sparsity) const
837  {
838  const unsigned int total_dofs = this->GetControlNDoFs();
839  dealii::CompressedSimpleSparsityPattern csp(total_dofs, total_dofs);
840 
841  dealii::DoFTools::make_sparsity_pattern(
842  static_cast<const dealii::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
843  csp);
844  this->GetControlDoFConstraints().condense(csp);
845  sparsity.copy_from(csp);
846  }
847 
851  template<>
852  void
853  DOpE::MethodOfLines_MultiMesh_SpaceTimeHandler<dealii::hp::FECollection,
854  dealii::hp::DoFHandler, dealii::BlockSparsityPattern,
855  dealii::BlockVector<double>, dope_dimension>::ComputeControlSparsityPattern(
856  dealii::BlockSparsityPattern & sparsity) const
857  {
858  const std::vector<unsigned int>& blocks = this->GetControlDoFsPerBlock();
859  dealii::BlockCompressedSimpleSparsityPattern csp(blocks.size(),
860  blocks.size());
861  for (unsigned int i = 0; i < blocks.size(); i++)
862  {
863  for (unsigned int j = 0; j < blocks.size(); j++)
864  {
865  csp.block(i, j).reinit(this->GetControlDoFsPerBlock(i),
866  this->GetControlDoFsPerBlock(j));
867  }
868  }
869  csp.collect_sizes();
870 
871  dealii::DoFTools::make_sparsity_pattern(
872  static_cast<const dealii::hp::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
873  csp);
874  this->GetControlDoFConstraints().condense(csp);
875  sparsity.copy_from(csp);
876  }
877 
878  /******************************************************/
879 
880  template<>
881  void
882  MethodOfLines_MultiMesh_SpaceTimeHandler<dealii::hp::FECollection,
883  dealii::hp::DoFHandler, dealii::SparsityPattern, dealii::Vector<double>,
884  dope_dimension>::ComputeControlSparsityPattern(
885  dealii::SparsityPattern & sparsity) const
886  {
887  const unsigned int total_dofs = this->GetControlNDoFs();
888  dealii::CompressedSimpleSparsityPattern csp(total_dofs, total_dofs);
889 
890  dealii::DoFTools::make_sparsity_pattern(
891  static_cast<const dealii::hp::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
892  csp);
893  this->GetControlDoFConstraints().condense(csp);
894  sparsity.copy_from(csp);
895  }
896 #endif //Endof explicit instanciation
897  /*******************************************************/
898  template<template<int, int> class FE, template<int, int> class DH,
899  typename SPARSITYPATTERN, typename VECTOR, int dim>
900  void
901  MethodOfLines_MultiMesh_SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
902  dim>::FlagIfLeftIsNotFinest(dealii::Triangulation<dim>& left,
903  const dealii::Triangulation<dim>& right)
904  {
905  auto element_list = GridTools::get_finest_common_cells(left, right);
906  auto element_iter = element_list.begin();
907  for (; element_iter != element_list.end(); element_iter++)
908  {
909  if (element_iter->second->has_children())
910  {
911  //left is not finest, so we should flag the left element to be refined
912  element_iter->first->set_refine_flag();
913  }
914  }
915  }
916 }
917 
918 #endif
void ReInitTime()
Definition: spacetimehandler_base.h:82
void SetActiveFEIndicesState(DOpEWrapper::DoFHandler< dealdim, DH > &dof_handler)
Definition: spacetimehandler.h:327
Definition: mapping_wrapper.h:48
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_multimesh_spacetimehandler.h:164
void RegisterMapping(const typename DOpEWrapper::Mapping< dealdim, DH > &mapping)
Definition: userdefineddofconstraints.h:76
Definition: dopetypes.h:50
void SetSparsityMaker(SparsityMaker< DH, dim > &sparsity_maker)
Definition: mol_multimesh_spacetimehandler.h:707
void RefineControlSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:481
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:500
void RefineStateSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:462
unsigned int GetConstraintNDoFs(std::string name) const
Definition: mol_multimesh_spacetimehandler.h:359
const DOpEWrapper::DoFHandler< dim, DH > & GetStateDoFHandler() const
Definition: mol_multimesh_spacetimehandler.h:233
void SetActiveFEIndicesControl(DOpEWrapper::DoFHandler< dopedim, DH > &dof_handler)
Definition: spacetimehandler.h:348
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, const dealii::Triangulation< 1 > &times, const Constraints &c, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dim, dim > &index_setter=ActiveFEIndexSetterInterface< dim, dim >())
Definition: mol_multimesh_spacetimehandler.h:119
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:284
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, const dealii::Triangulation< 1 > &times, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dim, dim > &index_setter=ActiveFEIndexSetterInterface< dim, dim >())
Definition: mol_multimesh_spacetimehandler.h:80
ControlType
Definition: dopetypes.h:102
void SetUserDefinedDoFConstraints(UserDefinedDoFConstraints< DH, dim > &constraints_maker)
Definition: mol_multimesh_spacetimehandler.h:692
Definition: dopetypes.h:50
void RefineSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:442
Definition: dopetypes.h:50
const std::vector< Point< dim > > & GetMapDoFToSupportPoints()
Definition: mol_multimesh_spacetimehandler.h:386
void ComputeControlSparsityPattern(SPARSITYPATTERN &sparsity) const
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dim, dim > &index_setter=ActiveFEIndexSetterInterface< dim, dim >())
Definition: mol_multimesh_spacetimehandler.h:61
Definition: timeiterator.h:63
unsigned int NewTimePointToOldTimePoint(unsigned int t) const
Definition: mol_multimesh_spacetimehandler.h:651
double get_k() const
Definition: timeiterator.h:195
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:95
RefinementType
Definition: dopetypes.h:48
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_multimesh_spacetimehandler.h:308
unsigned int GetControlNDoFs() const
Definition: mol_multimesh_spacetimehandler.h:343
const DOpEWrapper::DoFHandler< dim, DH > & GetControlDoFHandler() const
Definition: mol_multimesh_spacetimehandler.h:224
Definition: refinementcontainer.h:42
void ResetStateTriangulation(const dealii::Triangulation< dim > &tria)
Definition: mol_multimesh_spacetimehandler.h:722
const std::vector< unsigned int > & GetDoFsPerBlock(std::string name) const
Definition: constraints.h:159
unsigned int GetNGlobalConstraints() const
Definition: mol_multimesh_spacetimehandler.h:367
Definition: active_fe_index_setter_interface.h:39
unsigned int GetStateDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_multimesh_spacetimehandler.h:260
virtual const dealii::Vector< float > & GetLocalErrorIndicators() const
Definition: refinementcontainer.cc:46
virtual double GetTopFraction() const
Definition: refinementcontainer.cc:56
virtual void InterpolateState(VECTOR &result, const std::vector< VECTOR * > &local_vectors, double t, const TimeIterator &it) const
Definition: mol_multimesh_spacetimehandler.h:318
Definition: dopetypes.h:50
unsigned int GetConstraintDoFsPerBlock(std::string name, unsigned int b) const
Definition: mol_multimesh_spacetimehandler.h:268
void ComputeStateSparsityPattern(SPARSITYPATTERN &sparsity) const
Definition: mol_multimesh_spacetimehandler.h:400
Definition: constraints.h:33
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_multimesh_spacetimehandler.h:676
const std::vector< unsigned int > & GetControlDoFsPerBlock() const
Definition: mol_multimesh_spacetimehandler.h:276
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:85
const dealii::ConstraintMatrix & GetControlDoFConstraints() const
Definition: mol_multimesh_spacetimehandler.h:300
double get_right() const
Definition: timeiterator.h:185
Definition: spacetimehandler.h:71
void MapDoFsToSupportPoints(const DOpEWrapper::Mapping< dealdim, dealii::DoFHandler > &mapping, const DOpEWrapper::DoFHandler< dealdim, dealii::DoFHandler > &dh, VECTOR &support_points)
Definition: sth_internals.h:47
virtual void MakeControlDoFConstraints(const DOpEWrapper::DoFHandler< dopedim, DH > &dof_handler, dealii::ConstraintMatrix &dof_constraints) const
Definition: userdefineddofconstraints.h:103
const std::vector< unsigned int > & GetConstraintDoFsPerBlock(std::string name) const
Definition: mol_multimesh_spacetimehandler.h:292
void ReInit(std::vector< unsigned int > &control_dofs_per_block)
Definition: constraints.h:108
~MethodOfLines_MultiMesh_SpaceTimeHandler()
Definition: mol_multimesh_spacetimehandler.h:140
unsigned int GetControlDoFsPerBlock(unsigned int b, int=-1) const
Definition: mol_multimesh_spacetimehandler.h:252
Definition: mol_multimesh_spacetimehandler.h:57
Definition: sparsitymaker.h:50
bool UsesCoarsening() const
Definition: refinementcontainer.cc:103
double get_left() const
Definition: timeiterator.h:175
const FE< dim, dim > & GetFESystem(std::string name) const
Definition: mol_multimesh_spacetimehandler.h:413
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:351
void IncrementControlTicket()
Definition: spacetimehandler_base.h:382
Definition: solutiontransfer_wrapper.h:51
unsigned int GetNLocalConstraints() const
Definition: mol_multimesh_spacetimehandler.h:376
void RefineStateSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:515
Definition: dopeexception.h:35
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, const Constraints &c, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dim, dim > &index_setter=ActiveFEIndexSetterInterface< dim, dim >())
Definition: mol_multimesh_spacetimehandler.h:99
Definition: dopetypes.h:50
const DOpEWrapper::Mapping< dim, DH > & GetMapping() const
Definition: mol_multimesh_spacetimehandler.h:243
void SpatialMeshTransferControl(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_multimesh_spacetimehandler.h:664
void RefineControlSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:583
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
unsigned int n_dofs(std::string name) const
Definition: constraints.h:141
void IncrementStateTicket()
Definition: spacetimehandler_base.h:373
void ResetControlTriangulation(const dealii::Triangulation< dim > &tria)
Definition: mol_multimesh_spacetimehandler.h:742