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  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  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>&
276  GetControlDoFsPerBlock(int /*time_point*/= -1) const
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  InterpolateControl(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::InterpolateControl");
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 
339  virtual void
340  InterpolateState(VECTOR& result,
341  const std::vector<VECTOR*> & local_vectors, double t,
342  const TimeIterator& it) const
343  {
344  assert(it.get_left() <= t);
345  assert(it.get_right() >= t);
346  if (local_vectors.size() != 2)
347  throw DOpEException(
348  "This function is currently not implemented for anything other than"
349  " linear interpolation of 2 DoFs.",
350  "MethodOfLine_SpaceTimeHandler::InterpolateState");
351 
352  double lambda_l = (it.get_right() - t) / it.get_k();
353  double lambda_r = (t - it.get_left()) / it.get_k();
354 
355  //Here we assume that the numbering of dofs goes from left to right!
356  result = *local_vectors[0];
357 
358  result.sadd(lambda_l, lambda_r, *local_vectors[1]);
359  }
360 
364  unsigned int
365  GetControlNDoFs(int /*time_point*/= -1) const
366  {
367  return GetControlDoFHandler().n_dofs();
368  }
372  unsigned int
373  GetStateNDoFs(int /*time_point*/= -1) const
374  {
375  return GetStateDoFHandler().n_dofs();
376  }
380  unsigned int
381  GetConstraintNDoFs(std::string name) const
382  {
383  return constraints_.n_dofs(name);
384  }
388  unsigned int
390  {
391  return constraints_.n_dofs("global");
392  //return constraints_.global_n_dofs();
393  }
397  unsigned int
399  {
400  return constraints_.n_dofs("local");
401  //return constraints_.local_n_dofs();
402  }
403 
407  const std::vector<Point<dim> >&
409  {
410  support_points_.resize(GetStateNDoFs());
412  GetStateDoFHandler(), support_points_);
413  return support_points_;
414  }
415 
416  /******************************************************/
417  void
418  ComputeControlSparsityPattern(SPARSITYPATTERN & sparsity) const;
419 
420  /******************************************************/
421  void
422  ComputeStateSparsityPattern(SPARSITYPATTERN & sparsity) const
423  {
424  this->GetSparsityMaker()->ComputeSparsityPattern(
425  this->GetStateDoFHandler(), sparsity,
427  }
428 
429  /******************************************************/
430 
434  const FE<dim, dim>&
435  GetFESystem(std::string name) const
436  {
437  if (name == "state")
438  {
439  return *state_fe_;
440  }
441  else if (name == "control")
442  {
443  return *control_fe_;
444  }
445  else
446  {
447  throw DOpEException("Not implemented for name =" + name,
448  "MethodOfLines_MultiMesh_SpaceTimeHandler::GetFESystem");
449  }
450 
451  }
452 
463  void
466  {
467  assert(ref_type == DOpEtypes::RefinementType::global);
468  RefinementContainer ref_con_dummy;
469  RefineStateSpace(ref_con_dummy);
470  RefineControlSpace(ref_con_dummy);
471  }
472 
483  void
486  {
487  assert(ref_type == DOpEtypes::RefinementType::global);
488  RefinementContainer ref_con_dummy;
489  RefineStateSpace(ref_con_dummy);
490  }
491 
502  void
505  {
506  assert(ref_type == DOpEtypes::RefinementType::global);
507  RefinementContainer ref_con_dummy;
508  RefineControlSpace(ref_con_dummy);
509  }
510 
520  template<typename NUMBER>
521  void
522  RefineSpace(const RefinementContainer&ref_container)
523  {
524  RefineStateSpace(ref_container);
525  RefineControlSpace(ref_container);
526  }
527 
536  void
537  RefineStateSpace(const RefinementContainer& ref_container)
538  {
539  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
540 
541  //make sure that we do not use any coarsening
542  assert(!ref_container.UsesCoarsening());
543 
544  if (state_mesh_transfer_ != NULL)
545  {
546  delete state_mesh_transfer_;
547  state_mesh_transfer_ = NULL;
548  }
549  state_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dim, VECTOR,
550  DH>(state_dof_handler_);
551 
552  if (DOpEtypes::RefinementType::global == ref_type)
553  {
554  state_triangulation_.set_all_refine_flags();
555  }
556  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
557  {
558 
559  GridRefinement::refine_and_coarsen_fixed_number(
560  state_triangulation_, ref_container.GetLocalErrorIndicators(),
561  ref_container.GetTopFraction(),
562  ref_container.GetBottomFraction());
563  }
564  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
565  {
566 
567  GridRefinement::refine_and_coarsen_fixed_fraction(
568  state_triangulation_, ref_container.GetLocalErrorIndicators(),
569  ref_container.GetTopFraction(),
570  ref_container.GetBottomFraction());
571  }
572  else if (DOpEtypes::RefinementType::optimized == ref_type)
573  {
574 
575  GridRefinement::refine_and_coarsen_optimize(state_triangulation_,
576  ref_container.GetLocalErrorIndicators(),
577  ref_container.GetConvergenceOrder());
578  }
579  else if (DOpEtypes::RefinementType::finest_of_both == ref_type)
580  {
581  this->FlagIfLeftIsNotFinest(state_triangulation_,
582  control_triangulation_);
583  }
584  else
585  {
586  throw DOpEException("Not implemented for name =" + ref_type,
587  "MethodOfLines_MultiMesh_SpaceTimeHandler::RefineSpace");
588  }
589  state_triangulation_.prepare_coarsening_and_refinement();
590  if (state_mesh_transfer_ != NULL)
591  state_mesh_transfer_->prepare_for_pure_refinement();
592 
593  state_triangulation_.execute_coarsening_and_refinement();
594  }
595 
604  void
606  {
607  DOpEtypes::RefinementType ref_type = ref_container.GetRefType();
608 
609  //make sure that we do not use any coarsening
610  assert(!ref_container.UsesCoarsening());
611 
612  if (control_mesh_transfer_ != NULL)
613  {
614  delete control_mesh_transfer_;
615  control_mesh_transfer_ = NULL;
616  }
617 #if dope_dimension == deal_II_dimension
618  control_mesh_transfer_ = new DOpEWrapper::SolutionTransfer<dim, VECTOR,
619  DH>(control_dof_handler_);
620 #endif
621  if (DOpEtypes::RefinementType::global == ref_type)
622  {
623  control_triangulation_.set_all_refine_flags();
624  }
625 
626  else if (DOpEtypes::RefinementType::fixed_number == ref_type)
627  {
628 
629  GridRefinement::refine_and_coarsen_fixed_number(
630  control_triangulation_, ref_container.GetLocalErrorIndicators(),
631  ref_container.GetTopFraction(),
632  ref_container.GetBottomFraction());
633  }
634  else if (DOpEtypes::RefinementType::fixed_fraction == ref_type)
635  {
636 
637  GridRefinement::refine_and_coarsen_fixed_fraction(
638  control_triangulation_, ref_container.GetLocalErrorIndicators(),
639  ref_container.GetTopFraction(),
640  ref_container.GetBottomFraction());
641  }
642  else if (DOpEtypes::RefinementType::optimized == ref_type)
643  {
644 
645  GridRefinement::refine_and_coarsen_optimize(control_triangulation_,
646  ref_container.GetLocalErrorIndicators(),
647  ref_container.GetConvergenceOrder());
648  }
649  else if ("finest-of-both")
650  {
651  this->FlagIfLeftIsNotFinest(control_triangulation_,
652  state_triangulation_);
653  }
654  else
655  {
656  throw DOpEException("Not implemented for name =" + ref_type,
657  "MethodOfLines_MultiMesh_SpaceTimeHandler::RefineSpace");
658  }
659  control_triangulation_.prepare_coarsening_and_refinement();
660 
661  //FIXME: works only if no coarsening happens, because we do not have the vectors to be interpolated availiable...
662  if (control_mesh_transfer_ != NULL)
663  control_mesh_transfer_->prepare_for_pure_refinement();
664 
665  control_triangulation_.execute_coarsening_and_refinement();
666  }
667  /******************************************************/
668 
672  unsigned int
673  NewTimePointToOldTimePoint(unsigned int t) const
674  {
675  //TODO this has to be implemented when temporal refinement is possible!
676  //At present the temporal grid can't be refined
677  return t;
678  }
679 
680  /******************************************************/
681 
685  void
686  SpatialMeshTransferControl(const VECTOR& old_values,
687  VECTOR& new_values) const
688  {
689  if (control_mesh_transfer_ != NULL)
690  control_mesh_transfer_->refine_interpolate(old_values, new_values);
691  }
692  /******************************************************/
693 
697  void
698  SpatialMeshTransferState(const VECTOR& old_values,
699  VECTOR& new_values) const
700  {
701  if (state_mesh_transfer_ != NULL)
702  state_mesh_transfer_->refine_interpolate(old_values, new_values);
703  }
704  /******************************************************/
713  void
715  UserDefinedDoFConstraints<DH, dim>& constraints_maker)
716  {
717  user_defined_dof_constr_ = &constraints_maker;
718  user_defined_dof_constr_->RegisterMapping(this->GetMapping());
719  }
720  /******************************************************/
728  void
730  {
731  if (sparsitymaker_ != NULL && sparse_mkr_dynamic_)
732  delete sparsitymaker_;
733  sparsitymaker_ = &sparsity_maker;
734  sparse_mkr_dynamic_ = false;
735  }
736  /******************************************************/
743  void
744  ResetStateTriangulation(const dealii::Triangulation<dim>& tria)
745  {
746  state_dof_handler_.clear();
747  state_triangulation_.clear();
748  state_triangulation_.copy_triangulation(tria);
749  state_dof_handler_.initialize(state_triangulation_, *state_fe_);
750  this->IncrementStateTicket();
751  if (state_mesh_transfer_ != NULL)
752  delete state_mesh_transfer_;
753  state_mesh_transfer_ = NULL;
754  }
755 
756  /******************************************************/
763  void
764  ResetControlTriangulation(const dealii::Triangulation<dim>& tria)
765  {
766  control_dof_handler_.clear();
767  control_triangulation_.clear();
768  control_triangulation_.copy_triangulation(tria);
769  control_dof_handler_.initialize(control_triangulation_, *control_fe_);
770  this->IncrementControlTicket();
771  if (control_mesh_transfer_ != NULL)
772  delete control_mesh_transfer_;
773  control_mesh_transfer_ = NULL;
774  }
775 
776  private:
778  GetSparsityMaker() const
779  {
780  return sparsitymaker_;
781  }
782  const UserDefinedDoFConstraints<DH, dim>*
783  GetUserDefinedDoFConstraints() const
784  {
785  return user_defined_dof_constr_;
786  }
787 
788  void
789  FlagIfLeftIsNotFinest(dealii::Triangulation<dim>& left,
790  const dealii::Triangulation<dim>& right);
791 
792  SparsityMaker<DH, dim>* sparsitymaker_;
793  UserDefinedDoFConstraints<DH, dim>* user_defined_dof_constr_;
794 
795  dealii::Triangulation<dim> control_triangulation_;
796  dealii::Triangulation<dim>& state_triangulation_;
797  DOpEWrapper::DoFHandler<dim, DH> control_dof_handler_;
798  DOpEWrapper::DoFHandler<dim, DH> state_dof_handler_;
799 
800  std::vector<unsigned int> control_dofs_per_block_;
801  std::vector<unsigned int> state_dofs_per_block_;
802 
803  dealii::ConstraintMatrix control_dof_constraints_;
804  dealii::ConstraintMatrix state_dof_constraints_;
805 
806  const dealii::SmartPointer<const FE<dim, dim>> control_fe_;
807  const dealii::SmartPointer<const FE<dim, dim>> state_fe_;
808 
810 
811  std::vector<Point<dim> > support_points_;
812 
813  Constraints constraints_;
814  DOpEWrapper::SolutionTransfer<dim, VECTOR,DH>* control_mesh_transfer_;
815  DOpEWrapper::SolutionTransfer<dim, VECTOR,DH>* state_mesh_transfer_;
816  bool sparse_mkr_dynamic_;
817  };
818 
819  /**************************explicit instantiation*************/
820 #if dope_dimension == deal_II_dimension
821 
824  template<>
825  void
827  dealii::DoFHandler, dealii::BlockSparsityPattern,
828  dealii::BlockVector<double>, dope_dimension>::ComputeControlSparsityPattern(
829  dealii::BlockSparsityPattern & sparsity) const
830  {
831  const std::vector<unsigned int>& blocks = this->GetControlDoFsPerBlock();
832  dealii::BlockCompressedSimpleSparsityPattern csp(blocks.size(),
833  blocks.size());
834  for (unsigned int i = 0; i < blocks.size(); i++)
835  {
836  for (unsigned int j = 0; j < blocks.size(); j++)
837  {
838  csp.block(i, j).reinit(this->GetControlDoFsPerBlock(i),
839  this->GetControlDoFsPerBlock(j));
840  }
841  }
842  csp.collect_sizes();
843 
844  dealii::DoFTools::make_sparsity_pattern(
845  static_cast<const dealii::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
846  csp);
847  this->GetControlDoFConstraints().condense(csp);
848  sparsity.copy_from(csp);
849  }
850 
851  /******************************************************/
852 
853  template<>
854  void
856  dealii::DoFHandler, dealii::SparsityPattern, dealii::Vector<double>,
857  dope_dimension>::ComputeControlSparsityPattern(
858  dealii::SparsityPattern & sparsity) const
859  {
860  const unsigned int total_dofs = this->GetControlNDoFs();
861  dealii::CompressedSimpleSparsityPattern csp(total_dofs, total_dofs);
862 
863  dealii::DoFTools::make_sparsity_pattern(
864  static_cast<const dealii::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
865  csp);
866  this->GetControlDoFConstraints().condense(csp);
867  sparsity.copy_from(csp);
868  }
869 
873  template<>
874  void
875  DOpE::MethodOfLines_MultiMesh_SpaceTimeHandler<dealii::hp::FECollection,
876  dealii::hp::DoFHandler, dealii::BlockSparsityPattern,
877  dealii::BlockVector<double>, dope_dimension>::ComputeControlSparsityPattern(
878  dealii::BlockSparsityPattern & sparsity) const
879  {
880  const std::vector<unsigned int>& blocks = this->GetControlDoFsPerBlock();
881  dealii::BlockCompressedSimpleSparsityPattern csp(blocks.size(),
882  blocks.size());
883  for (unsigned int i = 0; i < blocks.size(); i++)
884  {
885  for (unsigned int j = 0; j < blocks.size(); j++)
886  {
887  csp.block(i, j).reinit(this->GetControlDoFsPerBlock(i),
888  this->GetControlDoFsPerBlock(j));
889  }
890  }
891  csp.collect_sizes();
892 
893  dealii::DoFTools::make_sparsity_pattern(
894  static_cast<const dealii::hp::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
895  csp);
896  this->GetControlDoFConstraints().condense(csp);
897  sparsity.copy_from(csp);
898  }
899 
900  /******************************************************/
901 
902  template<>
903  void
904  MethodOfLines_MultiMesh_SpaceTimeHandler<dealii::hp::FECollection,
905  dealii::hp::DoFHandler, dealii::SparsityPattern, dealii::Vector<double>,
906  dope_dimension>::ComputeControlSparsityPattern(
907  dealii::SparsityPattern & sparsity) const
908  {
909  const unsigned int total_dofs = this->GetControlNDoFs();
910  dealii::CompressedSimpleSparsityPattern csp(total_dofs, total_dofs);
911 
912  dealii::DoFTools::make_sparsity_pattern(
913  static_cast<const dealii::hp::DoFHandler<deal_II_dimension>&>(this->GetControlDoFHandler()),
914  csp);
915  this->GetControlDoFConstraints().condense(csp);
916  sparsity.copy_from(csp);
917  }
918 #endif //Endof explicit instanciation
919  /*******************************************************/
920  template<template<int, int> class FE, template<int, int> class DH,
921  typename SPARSITYPATTERN, typename VECTOR, int dim>
922  void
923  MethodOfLines_MultiMesh_SpaceTimeHandler<FE, DH, SPARSITYPATTERN, VECTOR,
924  dim>::FlagIfLeftIsNotFinest(dealii::Triangulation<dim>& left,
925  const dealii::Triangulation<dim>& right)
926  {
927  auto element_list = GridTools::get_finest_common_cells(left, right);
928  auto element_iter = element_list.begin();
929  for (; element_iter != element_list.end(); element_iter++)
930  {
931  if (element_iter->second->has_children())
932  {
933  //left is not finest, so we should flag the left element to be refined
934  element_iter->first->set_refine_flag();
935  }
936  }
937  }
938 }
939 
940 #endif
void ReInitTime()
Definition: spacetimehandler_base.h:86
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:52
unsigned int GetControlNDoFs(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:365
void SetSparsityMaker(SparsityMaker< DH, dim > &sparsity_maker)
Definition: mol_multimesh_spacetimehandler.h:729
void RefineControlSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:503
void RefineSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:522
void RefineStateSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:484
unsigned int GetConstraintNDoFs(std::string name) const
Definition: mol_multimesh_spacetimehandler.h:381
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
const std::vector< unsigned int > & GetStateDoFsPerBlock(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:284
ControlType
Definition: dopetypes.h:103
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, 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
void SetUserDefinedDoFConstraints(UserDefinedDoFConstraints< DH, dim > &constraints_maker)
Definition: mol_multimesh_spacetimehandler.h:714
Definition: dopetypes.h:52
void RefineSpace(DOpEtypes::RefinementType ref_type=DOpEtypes::RefinementType::global)
Definition: mol_multimesh_spacetimehandler.h:464
Definition: dopetypes.h:52
const std::vector< Point< dim > > & GetMapDoFToSupportPoints()
Definition: mol_multimesh_spacetimehandler.h:408
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:673
double get_k() const
Definition: timeiterator.h:195
DOpEtypes::RefinementType GetRefType() const
Definition: refinementcontainer.cc:95
RefinementType
Definition: dopetypes.h:50
const dealii::ConstraintMatrix & GetStateDoFConstraints() const
Definition: mol_multimesh_spacetimehandler.h:308
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:744
const std::vector< unsigned int > & GetDoFsPerBlock(std::string name) const
Definition: constraints.h:159
unsigned int GetNGlobalConstraints() const
Definition: mol_multimesh_spacetimehandler.h:389
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:340
Definition: dopetypes.h:52
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:422
Definition: constraints.h:33
void SpatialMeshTransferState(const VECTOR &old_values, VECTOR &new_values) const
Definition: mol_multimesh_spacetimehandler.h:698
virtual double GetConvergenceOrder() const
Definition: refinementcontainer.cc:85
const dealii::ConstraintMatrix & GetControlDoFConstraints() const
Definition: mol_multimesh_spacetimehandler.h:300
const std::vector< unsigned int > & GetControlDoFsPerBlock(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:276
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 InterpolateControl(VECTOR &result, const std::vector< VECTOR * > &local_vectors, double t, const TimeIterator &it) const
Definition: mol_multimesh_spacetimehandler.h:318
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:435
unsigned int GetStateNDoFs(int=-1) const
Definition: mol_multimesh_spacetimehandler.h:373
void IncrementControlTicket()
Definition: spacetimehandler_base.h:446
Definition: solutiontransfer_wrapper.h:51
unsigned int GetNLocalConstraints() const
Definition: mol_multimesh_spacetimehandler.h:398
void RefineStateSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:537
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
MethodOfLines_MultiMesh_SpaceTimeHandler(dealii::Triangulation< dim > &triangulation, const FE< dim, dim > &control_fe, const FE< dim, dim > &state_fe, dealii::Triangulation< 1 > &times, DOpEtypes::ControlType type, const ActiveFEIndexSetterInterface< dim, dim > &index_setter=ActiveFEIndexSetterInterface< dim, dim >())
Definition: mol_multimesh_spacetimehandler.h:80
Definition: dopetypes.h:52
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:686
void RefineControlSpace(const RefinementContainer &ref_container)
Definition: mol_multimesh_spacetimehandler.h:605
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:437
void ResetControlTriangulation(const dealii::Triangulation< dim > &tria)
Definition: mol_multimesh_spacetimehandler.h:764