DOpE
integrator_multimesh.h
Go to the documentation of this file.
1 
24 #ifndef IntegratorMultiMesh_H_
25 #define IntegratorMultiMesh_H_
26 
27 #include <lac/vector.h>
28 #include <lac/block_sparsity_pattern.h>
29 #include <lac/block_sparse_matrix.h>
30 
31 #include <numerics/matrix_tools.h>
32 #include <numerics/vector_tools.h>
33 #include <deal.II/grid/grid_tools.h>
34 #include <base/function.h>
35 
36 #include <vector>
37 
40 
41 namespace DOpE
42 {
60  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
61  int dim>
63  {
64  public:
65  IntegratorMultiMesh(INTEGRATORDATACONT& idc);
66 
68 
69  void
70  ReInit();
71 
72  template<typename PROBLEM>
73  void
74  ComputeNonlinearResidual(PROBLEM& pde, VECTOR &residual,
75  bool apply_boundary_values = true);
76  template<typename PROBLEM>
77  void
78  ComputeNonlinearLhs(PROBLEM& pde, VECTOR &residual,
79  bool apply_boundary_values = true);
80  template<typename PROBLEM>
81  void
82  ComputeNonlinearRhs(PROBLEM& pde, VECTOR &residual,
83  bool apply_boundary_values = true);
84  template<typename PROBLEM, typename MATRIX>
85  void
86  ComputeMatrix(PROBLEM& pde, MATRIX &matrix);
87 
88  template<typename PROBLEM>
89  void ComputeNonlinearAlgebraicResidual (PROBLEM& pde, VECTOR &residual);
90  template<typename PROBLEM>
91  void ComputeLocalControlConstraints (PROBLEM& pde, VECTOR &constraints);
92 
93  template<typename PROBLEM>
94  SCALAR
95  ComputeDomainScalar(PROBLEM& pde);
96  template<typename PROBLEM>
97  SCALAR
98  ComputePointScalar(PROBLEM& pde);
99  template<typename PROBLEM>
100  SCALAR
101  ComputeBoundaryScalar(PROBLEM& pde);
102  template<typename PROBLEM>
103  SCALAR
104  ComputeFaceScalar(PROBLEM& pde);
105  template<typename PROBLEM>
106  SCALAR
107  ComputeAlgebraicScalar(PROBLEM& pde);
108 
109  template<typename PROBLEM>
110  void
111  ApplyInitialBoundaryValues(PROBLEM& pde, VECTOR &u);
112  template<typename PROBLEM>
113  void
114  ApplyNewtonBoundaryValues(PROBLEM& pde, VECTOR &u);
115  template<typename PROBLEM, typename MATRIX>
116  void
117  ApplyNewtonBoundaryValues(PROBLEM& pde, MATRIX &matrix, VECTOR &rhs,
118  VECTOR &sol);
119 
120  inline void
121  AddDomainData(std::string name, const VECTOR* new_data);
122  inline void
123  DeleteDomainData(std::string name);
124 
125  inline void
126  AddParamData(std::string name, const dealii::Vector<SCALAR>* new_data);
127  inline void
128  DeleteParamData(std::string name);
129 
130  protected:
131  inline const std::map<std::string, const dealii::Vector<SCALAR>*>&
132  GetParamData() const;
133 
134  inline INTEGRATORDATACONT&
136  inline const std::map<std::string, const VECTOR*>&
137  GetDomainData() const;
138  inline void AddPresetRightHandSide(double s, VECTOR &residual) const;
139 
140  private:
141  template<template<int, int> class DH>
142  void
143  InterpolateBoundaryValues(
144  const DOpEWrapper::DoFHandler<dim, DH>* dof_handler,
145  const unsigned int color, const dealii::Function<dim>& function,
146  std::map<unsigned int, SCALAR>& boundary_values,
147  const std::vector<bool>& comp_mask) const;
148 
153  template<typename PROBLEM, template<int, int> class DH>
154  inline void ComputeNonlinearResidual_Recursive(
155  PROBLEM& pde,
156  VECTOR& residual,
157  typename std::vector<typename DH<dim, dim>::cell_iterator>& element_iter,
158  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_iter,
159  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
162 
167  template<typename PROBLEM, template<int, int> class DH>
168  inline void ComputeNonlinearRhs_Recursive(
169  PROBLEM& pde,
170  VECTOR& residual,
171  typename std::vector<typename DH<dim, dim>::cell_iterator>& element_iter,
172  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_iter,
173  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
176 
181  template<typename PROBLEM, typename MATRIX, template<int, int> class DH>
182  inline void ComputeMatrix_Recursive(
183  PROBLEM& pde,
184  MATRIX& matrix,
185  typename std::vector<typename DH<dim, dim>::cell_iterator>& element_iter,
186  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_iter,
187  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
190 
195  template<typename PROBLEM, template<int, int> class DH>
196  inline SCALAR ComputeDomainScalar_Recursive(
197  PROBLEM& pde,
198  typename std::vector<typename DH<dim, dim>::cell_iterator>& element_iter,
199  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_iter,
200  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
202 
207  template<typename PROBLEM, template<int, int> class DH>
208  inline SCALAR ComputeBoundaryScalar_Recursive(
209  PROBLEM& pde,
210  typename std::vector<typename DH<dim, dim>::cell_iterator>& element_iter,
211  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_iter,
212  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
214 
215  INTEGRATORDATACONT & idc_;
216 
217  std::map<std::string, const VECTOR*> domain_data_;
218  std::map<std::string, const dealii::Vector<SCALAR>*> param_data_;
219  };
220 
221  /**********************************Implementation*******************************************/
222 
223  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
224  int dim>
226  INTEGRATORDATACONT& idc) :
227  idc_(idc)
228  {
229  }
230 
231  /**********************************Implementation*******************************************/
232 
233  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
234  int dim>
236  {
237 
238  }
239 
240  /*******************************************************************************************/
241 
242  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
243  int dim>
244  void
246  {
247 
248  }
249 
250  /*******************************************************************************************/
251 
252  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
253  int dim>
254  template<typename PROBLEM>
255  void
257  PROBLEM& pde, VECTOR &residual, bool apply_boundary_values)
258  {
259  residual = 0.;
260 
261  const auto& dof_handler =
262  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
263 
264  assert(dof_handler.size() == 2);
265 
266  const auto tria_element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler().get_tria(),
267  dof_handler[1]->GetDEALDoFHandler().get_tria());
268  const auto element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler(),
269  dof_handler[1]->GetDEALDoFHandler());
270 
271  auto element_iter = element_list.begin();
272  auto tria_element_iter = tria_element_list.begin();
273 
274  std::vector<decltype(tria_element_iter->first)> tria_element(2);
275  tria_element[0] = tria_element_iter->first;
276  tria_element[1] = tria_element_iter->second;
277 
278  std::vector<decltype(element_iter->first)> element(2);
279  element[0] = element_iter->first;
280  element[1] = element_iter->second;
281  int coarse_index = 0; //element[coarse_index] is the coarser of the two.
282  int fine_index = 0; //element[fine_index] is the finer of the two (both indices are 2 = element.size() if they are equally refined.
283 
284  // Generate the data containers.
285  idc_.InitializeMMEDC(pde.GetUpdateFlags(),
286  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
287  this->GetParamData(), this->GetDomainData());
288  auto& edc = idc_.GetMultimeshElementDataContainer();
289 
290  bool need_interfaces = pde.HasInterfaces();
291  idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
292  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
293  element,
294  tria_element,
295  this->GetParamData(),
296  this->GetDomainData(),
297  need_interfaces);
298  auto& fdc = idc_.GetMultimeshFaceDataContainer();
299 
300  for(; element_iter != element_list.end(); element_iter++)
301  {
302  element[0] = element_iter->first;
303  element[1] = element_iter->second;
304  tria_element[0] = tria_element_iter->first;
305  tria_element[1] = tria_element_iter->second;
306  FullMatrix<SCALAR> prolong_matrix;
307 
308  if(element[0]->has_children())
309  {
310  prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
311  coarse_index =1;
312  fine_index = 0;
313  }
314  else
315  {
316  if(element[1]->has_children())
317  {
318  prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
319  coarse_index =0;
320  fine_index = 1;
321  }
322  else
323  {
324  assert(element.size() ==2);
325  coarse_index = fine_index = 2;
326  }
327  }
328  ComputeNonlinearResidual_Recursive(pde,residual,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
329  tria_element_iter++;
330  }
331  //Check if some preset righthandside exists.
332  AddPresetRightHandSide(-1.,residual);
333 
334  if (apply_boundary_values)
335  {
336  ApplyNewtonBoundaryValues(pde, residual);
337  }
338  }
339 
340  /*******************************************************************************************/
341 
342  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
343  int dim>
344  template<typename PROBLEM>
345  void
347  PROBLEM& pde, VECTOR &residual, bool apply_boundary_values)
348  {
349  {
350  throw DOpEException("This function needs to be implemented!", "IntegratorMultiMesh::ComputeNonlinearLhs");
351  }
352  }
353 
354  /*******************************************************************************************/
355 
356  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
357  int dim>
358  template<typename PROBLEM>
359  void
361  PROBLEM& pde, VECTOR &residual, bool apply_boundary_values)
362  {
363  residual = 0.;
364 
365  const auto& dof_handler =
366  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
367 
368  assert(dof_handler.size() == 2);
369 
370  const auto tria_element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler().get_tria(),
371  dof_handler[1]->GetDEALDoFHandler().get_tria());
372  const auto element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler(),
373  dof_handler[1]->GetDEALDoFHandler());
374 
375  auto element_iter = element_list.begin();
376  auto tria_element_iter = tria_element_list.begin();
377 
378  std::vector<decltype(tria_element_iter->first)> tria_element(2);
379  tria_element[0] = tria_element_iter->first;
380  tria_element[1] = tria_element_iter->second;
381 
382  std::vector<decltype(element_iter->first)> element(2);
383  element[0] = element_iter->first;
384  element[1] = element_iter->second;
385  int coarse_index = 0; //element[coarse_index] is the coarser of the two.
386  int fine_index = 0; //element[fine_index] is the finer of the two (both indices are 2 = element.size() if they are equally refined.
387 
388  // Generate the data containers.
389  idc_.InitializeMMEDC(pde.GetUpdateFlags(),
390  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
391  this->GetParamData(), this->GetDomainData());
392  auto& edc = idc_.GetMultimeshElementDataContainer();
393 
394  idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
395  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
396  element,
397  tria_element,
398  this->GetParamData(),
399  this->GetDomainData());
400  auto& fdc = idc_.GetMultimeshFaceDataContainer();
401 
402  for(; element_iter != element_list.end(); element_iter++)
403  {
404  element[0] = element_iter->first;
405  element[1] = element_iter->second;
406  tria_element[0] = tria_element_iter->first;
407  tria_element[1] = tria_element_iter->second;
408  FullMatrix<SCALAR> prolong_matrix;
409 
410  if(element[0]->has_children())
411  {
412  prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
413  coarse_index =1;
414  fine_index = 0;
415  }
416  else
417  {
418  if(element[1]->has_children())
419  {
420  prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
421  coarse_index =0;
422  fine_index = 1;
423  }
424  else
425  {
426  assert(element.size() ==2);
427  coarse_index = fine_index = 2;
428  }
429  }
430  ComputeNonlinearRhs_Recursive(pde,residual,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
431  tria_element_iter++;
432  }
433  //Check if some preset righthandside exists.
434  AddPresetRightHandSide(1.,residual);
435 
436  if (apply_boundary_values)
437  {
438  ApplyNewtonBoundaryValues(pde, residual);
439  }
440  }
441 
442  /*******************************************************************************************/
443 
444  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
445  int dim>
446  template<typename PROBLEM, typename MATRIX>
447  void
449  PROBLEM& pde, MATRIX &matrix)
450  {
451  matrix = 0.;
452 
453  const auto& dof_handler =
454  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
455 
456  assert(dof_handler.size() == 2);
457 
458  const auto tria_element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler().get_tria(),
459  dof_handler[1]->GetDEALDoFHandler().get_tria());
460  const auto element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler(),
461  dof_handler[1]->GetDEALDoFHandler());
462 
463  auto element_iter = element_list.begin();
464  auto tria_element_iter = tria_element_list.begin();
465 
466  std::vector<decltype(tria_element_iter->first)> tria_element(2);
467  tria_element[0] = tria_element_iter->first;
468  tria_element[1] = tria_element_iter->second;
469 
470  std::vector<decltype(element_iter->first)> element(2);
471  element[0] = element_iter->first;
472  element[1] = element_iter->second;
473  int coarse_index = 0; //element[coarse_index] is the coarser of the two.
474  int fine_index = 0; //element[fine_index] is the finer of the two (both indices are 2 = element.size() if they are equally refined.
475 
476  // Generate the data containers.
477  idc_.InitializeMMEDC(pde.GetUpdateFlags(),
478  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
479  this->GetParamData(), this->GetDomainData());
480  auto& edc = idc_.GetMultimeshElementDataContainer();
481 
482  bool need_interfaces = pde.HasInterfaces();
483  idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
484  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
485  element,
486  tria_element,
487  this->GetParamData(),
488  this->GetDomainData(),
489  need_interfaces);
490  auto & fdc = idc_.GetMultimeshFaceDataContainer();
491 
492  for(; element_iter != element_list.end(); element_iter++)
493  {
494  element[0] = element_iter->first;
495  element[1] = element_iter->second;
496  tria_element[0] = tria_element_iter->first;
497  tria_element[1] = tria_element_iter->second;
498  FullMatrix<SCALAR> prolong_matrix;
499 
500  if(element[0]->has_children())
501  {
502  prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
503  coarse_index =1;
504  fine_index = 0;
505  }
506  else
507  {
508  if(element[1]->has_children())
509  {
510  prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
511  coarse_index =0;
512  fine_index = 1;
513  }
514  else
515  {
516  assert(element.size() ==2);
517  coarse_index = fine_index = 2;
518  }
519  }
520  ComputeMatrix_Recursive(pde,matrix,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
521  tria_element_iter++;
522  }
523  }
524 
525  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
526  int dim>
527  template<typename PROBLEM>
528  SCALAR
530  PROBLEM& pde)
531  {
532  SCALAR ret = 0.;
533 
534  const auto& dof_handler =
535  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
536 
537  assert(dof_handler.size() == 2);
538 
539  const auto tria_element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler().get_tria(),
540  dof_handler[1]->GetDEALDoFHandler().get_tria());
541  const auto element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler(),
542  dof_handler[1]->GetDEALDoFHandler());
543 
544  auto element_iter = element_list.begin();
545  auto tria_element_iter = tria_element_list.begin();
546 
547  std::vector<decltype(tria_element_iter->first)> tria_element(2);
548  tria_element[0] = tria_element_iter->first;
549  tria_element[1] = tria_element_iter->second;
550 
551  std::vector<decltype(element_iter->first)> element(2);
552  element[0] = element_iter->first;
553  element[1] = element_iter->second;
554  int coarse_index = 0; //element[coarse_index] is the coarser of the two.
555  int fine_index = 0; //element[fine_index] is the finer of the two (both indices are 2 = element.size() if they are equally refined.
556 
557  if (pde.HasFaces())
558  {
559  throw DOpEException("This function should not be called when faces are needed!",
560  "IntegratorMultiMesh::ComputeDomainScalar");
561  }
562 
563  // Generate the data containers.
564  idc_.InitializeMMEDC(pde.GetUpdateFlags(),
565  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
566  this->GetParamData(), this->GetDomainData());
567  auto& edc = idc_.GetMultimeshElementDataContainer();
568 
569  for(; element_iter != element_list.end(); element_iter++)
570  {
571  element[0] = element_iter->first;
572  element[1] = element_iter->second;
573  tria_element[0] = tria_element_iter->first;
574  tria_element[1] = tria_element_iter->second;
575  FullMatrix<SCALAR> prolong_matrix;
576 
577  if(element[0]->has_children())
578  {
579  prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
580  coarse_index =1;
581  fine_index = 0;
582  }
583  else
584  {
585  if(element[1]->has_children())
586  {
587  prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
588  coarse_index =0;
589  fine_index = 1;
590  }
591  else
592  {
593  assert(element.size() ==2);
594  coarse_index = fine_index = 2;
595  }
596  }
597  ret += ComputeDomainScalar_Recursive(pde,element,tria_element,prolong_matrix,coarse_index,fine_index,edc);
598  tria_element_iter++;
599  }
600  return ret;
601  }
602 
603 
604  /*******************************************************************************************/
605 
606  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
607  int dim>
608  template<typename PROBLEM>
609  SCALAR
611  PROBLEM& pde)
612  {
613 
614  {
615  SCALAR ret = 0.;
616  ret += pde.PointFunctional(this->GetParamData(),
617  this->GetDomainData());
618 
619  return ret;
620  }
621  }
622  /*******************************************************************************************/
623 
624  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
625  int dim>
626  template<typename PROBLEM>
627  SCALAR
629  PROBLEM& pde
630  )
631  {
632  SCALAR ret = 0.;
633  // Begin integration
634  const auto& dof_handler =
635  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
636 
637  const auto tria_element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler().get_tria(),
638  dof_handler[1]->GetDEALDoFHandler().get_tria());
639  const auto element_list = GridTools::get_finest_common_cells (dof_handler[0]->GetDEALDoFHandler(),
640  dof_handler[1]->GetDEALDoFHandler());
641 
642  auto element_iter = element_list.begin();
643  auto tria_element_iter = tria_element_list.begin();
644 
645  std::vector<decltype(tria_element_iter->first)> tria_element(2);
646  tria_element[0] = tria_element_iter->first;
647  tria_element[1] = tria_element_iter->second;
648 
649  std::vector<decltype(element_iter->first)> element(2);
650  element[0] = element_iter->first;
651  element[1] = element_iter->second;
652  int coarse_index = 0; //element[coarse_index] is the coarser of the two.
653  int fine_index = 0; //element[fine_index] is the finer of the two (both indices are 2 = element.size() if they are equally refined.
654 
655 
656  idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
657  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
658  element, tria_element,
659  this->GetParamData(),
660  this->GetDomainData());
661  auto & fdc = idc_.GetMultimeshFaceDataContainer();
662 
663  std::vector<unsigned int> boundary_functional_colors = pde.GetBoundaryFunctionalColors();
664  bool need_boundary_integrals = (boundary_functional_colors.size() > 0);
665  if(!need_boundary_integrals)
666  {
667  throw DOpEException("No boundary colors given!","IntegratorMultiMesh::ComputeBoundaryScalar");
668  }
669 
670  for(; element_iter != element_list.end(); element_iter++)
671  {
672  element[0] = element_iter->first;
673  element[1] = element_iter->second;
674  tria_element[0] = tria_element_iter->first;
675  tria_element[1] = tria_element_iter->second;
676  FullMatrix<SCALAR> prolong_matrix;
677 
678  if(element[0]->has_children())
679  {
680  prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
681  coarse_index =1;
682  fine_index = 0;
683  }
684  else
685  {
686  if(element[1]->has_children())
687  {
688  prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
689  coarse_index =0;
690  fine_index = 1;
691  }
692  else
693  {
694  assert(element.size() ==2);
695  coarse_index = fine_index = 2;
696  }
697  }
698  ret += ComputeBoundaryScalar_Recursive(pde,element,tria_element,prolong_matrix,coarse_index,fine_index,fdc);
699  tria_element_iter++;
700  }
701 
702  return ret;
703  }
704  /*******************************************************************************************/
705 
706  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
707  int dim>
708  template<typename PROBLEM>
709  SCALAR
711  PROBLEM& /*pde*/
712  )
713  {
714  throw DOpEException("This function needs to be implemented!", "IntegratorMultiMesh::ComputeFaceScalar");
715 // {
716 // SCALAR ret = 0.;
717 //#if deal_II_dimension == 2 || deal_II_dimension == 3
718 // // Begin integration
719 // const auto& dof_handler =
720 // pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
721 // auto element = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerBeginActive();
722 // auto endc = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerEnd();
723 //
724 // idc_.InitializeFDC(pde.GetFaceUpdateFlags(),
725 // *(pde.GetBaseProblem().GetSpaceTimeHandler()),
726 // element,
727 // this->GetParamData(),
728 // this->GetDomainData());
729 // auto & fdc = idc_.GetFaceDataContainer();
730 //
731 // bool need_faces = pde.HasFaces();
732 // if(!need_faces)
733 // {
734 // throw DOpEException("No faces required!","IntegratorMultiMesh::ComputeFaceScalar");
735 // }
736 //
737 // for (;element[0]!=endc[0]; element[0]++)
738 // {
739 // for(unsigned int dh=1; dh<dof_handler.size(); dh++)
740 // {
741 // if( element[dh] == endc[dh])
742 // {
743 // throw DOpEException("Elementnumbers in DoFHandlers are not matching!","IntegratorMultiMesh::ComputeFaceScalar");
744 // }
745 // }
746 //
747 // if(need_faces)
748 // {
749 // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
750 // {
751 // if (element[0]->neighbor_index(face) != -1)
752 // {
753 // fdc.ReInit(face);
754 // ret +=pde.FaceFunctional(fdc);
755 // }
756 // }
757 // }
758 // for(unsigned int dh=1; dh<dof_handler.size(); dh++)
759 // {
760 // element[dh]++;
761 // }
762 // }
763 //#else
764 // throw DOpEException("Not implemented in this dimension!",
765 // "IntegratorMultiMesh::ComputeFaceScalar");
766 //#endif
767 // return ret;
768 // }
769  }
770  /*******************************************************************************************/
771 
772  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
773  int dim>
774  template<typename PROBLEM>
775  SCALAR
777  PROBLEM& pde)
778  {
779 
780  {
781  SCALAR ret = 0.;
782  ret = pde.AlgebraicFunctional(this->GetParamData(),
783  this->GetDomainData());
784  return ret;
785  }
786  }
787 
788  /*******************************************************************************************/
789 
790  template <typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
791  template<typename PROBLEM>
793  ::ComputeNonlinearAlgebraicResidual (PROBLEM& pde, VECTOR &residual)
794  {
795  residual = 0.;
796  pde.AlgebraicResidual(residual,this->GetParamData(),this->GetDomainData());
797  //Check if some preset righthandside exists.
798  AddPresetRightHandSide(-1.,residual);
799  }
800 
801  /*******************************************************************************************/
802 
803  template <typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
804  template<typename PROBLEM>
806  ::ComputeLocalControlConstraints (PROBLEM& pde, VECTOR &constraints)
807  {
808  constraints = 0.;
809  pde.ComputeLocalControlConstraints(constraints,this->GetParamData(),this->GetDomainData());
810  }
811 
812  /*******************************************************************************************/
813 
814  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
815  int dim>
816  template<typename PROBLEM>
817  void
819  PROBLEM& pde, VECTOR &u)
820  {
821  //TODO Apply constraints locally, see, e.g., dealii step-27 ? But howto do this in the newton iter
822  // e.g. sometimes we need zero sometimes we need other values.
823 
824  //Never Condense Nodes Here ! Or All will fail if the state is not initialized with zero!
825  //pde.GetDoFConstraints().condense(u);
826  std::vector<unsigned int> dirichlet_colors = pde.GetDirichletColors();
827  for (unsigned int i = 0; i < dirichlet_colors.size(); i++)
828  {
829  unsigned int color = dirichlet_colors[i];
830  std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
831  std::map<unsigned int, SCALAR> boundary_values;
832 
833  InterpolateBoundaryValues(pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0], color, pde.GetDirichletValues(color, this->GetParamData(),
834  this->GetDomainData()), boundary_values, comp_mask );
835 
836  for (typename std::map<unsigned int, SCALAR>::const_iterator p =
837  boundary_values.begin(); p != boundary_values.end(); p++)
838  {
839  u(p->first) = p->second;
840  }
841  }
842  }
843 
844  /*******************************************************************************************/
845 
846  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
847  int dim>
848  template<typename PROBLEM>
849  void
851  PROBLEM& pde, VECTOR &u)
852  {
853  //TODO Apply constraints locally, see, e.g., dealii step-27 ? But howto do this in the newton iter
854  // e.g. sometimes we need zero sometimes we need other values.
855 
856  pde.GetDoFConstraints().condense(u);
857  std::vector<unsigned int> dirichlet_colors = pde.GetDirichletColors();
858  for (unsigned int i = 0; i < dirichlet_colors.size(); i++)
859  {
860  unsigned int color = dirichlet_colors[i];
861  std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
862  std::map<unsigned int, SCALAR> boundary_values;
863 
864 
865  InterpolateBoundaryValues(
866  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0],
867  color, dealii::ZeroFunction<dim>(comp_mask.size()),
868  boundary_values, comp_mask);
869 
870  for (typename std::map<unsigned int, SCALAR>::const_iterator p =
871  boundary_values.begin(); p != boundary_values.end(); p++)
872  {
873  u(p->first) = p->second;
874  }
875  }
876  }
877  /*******************************************************************************************/
878 
879  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
880  int dim>
881  template<typename PROBLEM, typename MATRIX>
882  void
884  PROBLEM& pde, MATRIX& matrix, VECTOR &rhs, VECTOR &sol)
885  {
886  //TODO Apply constraints locally, see, e.g., dealii step-27 ? But howto do this in the newton iter
887  // e.g. sometimes we need zero sometimes we need other values.
888  pde.GetDoFConstraints().condense(rhs);
889  pde.GetDoFConstraints().condense(matrix);
890  std::vector<unsigned int> dirichlet_colors = pde.GetDirichletColors();
891  for (unsigned int i = 0; i < dirichlet_colors.size(); i++)
892  {
893  unsigned int color = dirichlet_colors[i];
894  std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
895  std::map<unsigned int, SCALAR> boundary_values;
896 
897  InterpolateBoundaryValues(
898  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0],
899  color, dealii::ZeroFunction<dim>(comp_mask.size()),
900  boundary_values, comp_mask);
901 
902  dealii::MatrixTools::apply_boundary_values(boundary_values, matrix,
903  sol, rhs);
904  }
905  }
906 
907  /*******************************************************************************************/
908 
909  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
910  int dim>
911  void
913  std::string name, const VECTOR* new_data)
914  {
915  if (domain_data_.find(name) != domain_data_.end())
916  {
917  throw DOpEException(
918  "Adding multiple Data with name " + name + " is prohibited!",
919  "IntegratorMultiMesh::AddDomainData");
920  }
921  domain_data_.insert(std::pair<std::string, const VECTOR*>(name, new_data));
922  }
923 
924  /*******************************************************************************************/
925 
926  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
927  int dim>
928  void
930  std::string name)
931  {
932  typename std::map<std::string, const VECTOR *>::iterator it =
933  domain_data_.find(name);
934  if (it == domain_data_.end())
935  {
936  throw DOpEException(
937  "Deleting Data " + name + " is impossible! Data not found",
938  "IntegratorMultiMesh::DeleteDomainData");
939  }
940  domain_data_.erase(it);
941  }
942  /*******************************************************************************************/
943 
944  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
945  int dim>
946  void
948  VECTOR &residual) const
949  {
950  typename std::map<std::string, const VECTOR *>::const_iterator it =
951  domain_data_.find("fixed_rhs");
952  if (it != domain_data_.end())
953  {
954  assert(residual.size() == it->second->size());
955  residual.add(s,*(it->second));
956  }
957  }
958 
959  /*******************************************************************************************/
960 
961  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
962  int dim>
963  const std::map<std::string, const VECTOR*>&
965  {
966  return domain_data_;
967  }
968 
969  /*******************************************************************************************/
970 
971  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
972  int dim>
973  void
975  std::string name, const dealii::Vector<SCALAR>* new_data)
976  {
977  if (param_data_.find(name) != param_data_.end())
978  {
979  throw DOpEException(
980  "Adding multiple Data with name " + name + " is prohibited!",
981  "IntegratorMultiMesh::AddParamData");
982  }
983  param_data_.insert(
984  std::pair<std::string, const dealii::Vector<SCALAR>*>(name, new_data));
985  }
986 
987  /*******************************************************************************************/
988 
989  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
990  int dim>
991  void
993  std::string name)
994  {
995  typename std::map<std::string, const dealii::Vector<SCALAR>*>::iterator
996  it = param_data_.find(name);
997  if (it == param_data_.end())
998  {
999  throw DOpEException(
1000  "Deleting Data " + name + " is impossible! Data not found",
1001  "IntegratorMultiMesh::DeleteParamData");
1002  }
1003  param_data_.erase(it);
1004  }
1005 
1006  /*******************************************************************************************/
1007 
1008  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
1009  int dim>
1010  const std::map<std::string, const dealii::Vector<SCALAR>*>&
1012  {
1013  return param_data_;
1014  }
1015 
1016  /*******************************************************************************************/
1017 
1018  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
1019  int dim>
1020  INTEGRATORDATACONT&
1022  {
1023  return idc_;
1024  }
1025 
1026 
1027  /*******************************************************************************************/
1028 
1029  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
1030  int dim>
1031  template<template<int, int> class DH>
1032  void
1034  const DOpEWrapper::DoFHandler<dim, DH>* dof_handler,
1035  const unsigned int color, const dealii::Function<dim>& function,
1036  std::map<unsigned int, SCALAR>& boundary_values,
1037  const std::vector<bool>& comp_mask) const
1038  {
1039  dealii::VectorTools::interpolate_boundary_values(
1040  dof_handler->GetDEALDoFHandler(), color, function,
1041  boundary_values, comp_mask);
1042  }
1043 
1044  /*******************************************************************************************/
1045 
1046  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
1047  template<typename PROBLEM, template<int, int> class DH>
1048  void IntegratorMultiMesh<INTEGRATORDATACONT, VECTOR, SCALAR, dim>::ComputeNonlinearResidual_Recursive(
1049  PROBLEM& pde,
1050  VECTOR& residual,
1051  typename std::vector<typename DH<dim, dim>::cell_iterator> &element,
1052  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator> &tria_element,
1053  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
1054  Multimesh_ElementDataContainer<DH, VECTOR, dim>& edc,Multimesh_FaceDataContainer<DH, VECTOR, dim>& fdc)
1055  {
1056  if(!element[0]->has_children() && ! element[1]->has_children())
1057  {
1058  unsigned int dofs_per_element;
1059  dealii::Vector<SCALAR> local_vector;
1060  std::vector<unsigned int> local_dof_indices;
1061 
1062  bool need_faces = pde.HasFaces();
1063  std::vector<unsigned int> boundary_equation_colors = pde.GetBoundaryEquationColors();
1064  bool need_boundary_integrals = (boundary_equation_colors.size() > 0);
1065  bool need_interfaces = pde.HasInterfaces();
1066 
1067  edc.ReInit(coarse_index,fine_index,prolong_matrix);
1068 
1069  dofs_per_element = element[0]->get_fe().dofs_per_cell;
1070 
1071  local_vector.reinit(dofs_per_element);
1072  local_vector = 0;
1073 
1074  local_dof_indices.resize(0);
1075  local_dof_indices.resize(dofs_per_element, 0);
1076 
1077  //the second '1' plays only a role in the stationary case. In the non-stationary
1078  //case, scale_ico is set by the time-stepping-scheme
1079  pde.ElementEquation(edc, local_vector, 1., 1.);
1080  pde.ElementRhs(edc, local_vector, -1.);
1081 
1082  //FIXME Integrate on Faces of Element[0] that contain a fine-element face.
1083  if(need_faces || need_interfaces)
1084  {
1085  throw DOpEException(" Faces on multiple meshes not implemented yet!",
1086  "IntegratorMultiMesh::ComputeNonlinearResidual_Recursive");
1087  }
1088  unsigned int b_index = fine_index%2; //This takes care that if fine_index ==2 then we select the
1089  //zeros entry in the element vector
1090 
1091  if(need_boundary_integrals && element[b_index]->at_boundary())
1092  {
1093  for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1094  {
1095  if (element[b_index]->face(face)->at_boundary()
1096  &&
1097  (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1098  element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
1099  {
1100  fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1101  pde.BoundaryEquation(fdc,local_vector, 1., 1.);
1102  pde.BoundaryRhs(fdc,local_vector,-1.);
1103  }
1104  }
1105  }
1106  // if(need_faces)
1107  // {
1108  // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1109  // {
1110  // if (element[0]->neighbor_index(face) != -1)
1111  // {
1112  // fdc.ReInit(face);
1113  // pde.FaceEquation(fdc, local_vector);
1114  // pde.FaceRhs(fdc, local_vector,-1.);
1115  // }
1116  // }
1117  // }
1118  // if( need_interfaces)
1119  // {
1120  // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1121  // {
1122  // fdc.ReInit(face);
1123  // if (element[0]->neighbor_index(face) != -1
1124  // &&
1125  // fdc.GetMaterialId()!= fdc.GetNbrMaterialId())
1126  // {
1127  // fdc.ReInitNbr();
1128  // pde.InterfaceEquation(fdc, local_vector);
1129  // }
1130  // }
1131  // }
1132  if(coarse_index == 0) //Need to transfer the computed residual to the one on the coarse element
1133  {
1134  dealii::Vector<SCALAR> tmp(dofs_per_element);
1135  prolong_matrix.Tvmult(tmp,local_vector);
1136 
1137  //LocalToGlobal
1138  element[0]->get_dof_indices(local_dof_indices);
1139  for (unsigned int i = 0; i < dofs_per_element; ++i)
1140  {
1141  residual(local_dof_indices[i]) += tmp(i);
1142  }
1143  }
1144  else //Testfunctions are already the right ones...
1145  {
1146  //LocalToGlobal
1147  element[0]->get_dof_indices(local_dof_indices);
1148  for (unsigned int i = 0; i < dofs_per_element; ++i)
1149  {
1150  residual(local_dof_indices[i]) += local_vector(i);
1151  }
1152  }
1153 
1154  }//Endof the case on the finest level
1155  else
1156  {
1157  assert(fine_index != coarse_index);
1158  assert(element[fine_index]->has_children());
1159  assert(!element[coarse_index]->has_children());
1160  assert(tria_element[fine_index]->has_children());
1161  assert(!tria_element[coarse_index]->has_children());
1162 
1163  unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
1164 
1165  typename DH<dim, dim>::cell_iterator dofh_fine = element[fine_index];
1166  typename dealii::Triangulation<dim>::cell_iterator tria_fine = tria_element[fine_index];
1167 
1168  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1169  {
1170  FullMatrix<SCALAR> new_matrix(local_n_dofs);
1171  element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1172  prolong_matrix);
1173  element[fine_index] = dofh_fine->child(child);
1174  tria_element[fine_index] = tria_fine->child(child);
1175 
1176  ComputeNonlinearResidual_Recursive(pde,residual,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
1177  }
1178  }
1179  }
1180  /*******************************************************************************************/
1181 
1182  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
1183  template<typename PROBLEM, template<int, int> class DH>
1184  void IntegratorMultiMesh<INTEGRATORDATACONT, VECTOR, SCALAR, dim>::ComputeNonlinearRhs_Recursive(
1185  PROBLEM& pde,
1186  VECTOR& residual,
1187  typename std::vector<typename DH<dim, dim>::cell_iterator> &element,
1188  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator> &tria_element,
1189  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
1190  Multimesh_ElementDataContainer<DH, VECTOR, dim>& edc,Multimesh_FaceDataContainer<DH, VECTOR, dim>& fdc)
1191  {
1192  if(!element[0]->has_children() && ! element[1]->has_children())
1193  {
1194  unsigned int dofs_per_element;
1195  dealii::Vector<SCALAR> local_vector;
1196  std::vector<unsigned int> local_dof_indices;
1197 
1198  bool need_faces = pde.HasFaces();
1199  std::vector<unsigned int> boundary_equation_colors = pde.GetBoundaryEquationColors();
1200  bool need_boundary_integrals = (boundary_equation_colors.size() > 0);
1201  bool need_interfaces = pde.HasInterfaces();
1202 
1203  edc.ReInit(coarse_index,fine_index,prolong_matrix);
1204 
1205  dofs_per_element = element[0]->get_fe().dofs_per_cell;
1206 
1207  local_vector.reinit(dofs_per_element);
1208  local_vector = 0;
1209 
1210  local_dof_indices.resize(0);
1211  local_dof_indices.resize(dofs_per_element, 0);
1212 
1213  //the second '1' plays only a role in the stationary case. In the non-stationary
1214  //case, scale_ico is set by the time-stepping-scheme
1215  pde.ElementRhs(edc, local_vector, 1.);
1216 
1217  //FIXME Integrate on Faces of Element[0] that contain a fine-element face.
1218  if(need_faces )
1219  {
1220  throw DOpEException(" Faces on multiple meshes not implemented yet!",
1221  "IntegratorMultiMesh::ComputeNonlinearRhs_Recursive");
1222  }
1223  if(need_interfaces )
1224  {
1225  throw DOpEException(" Faces on multiple meshes not implemented yet!",
1226  "IntegratorMultiMesh::ComputeNonlinearRhs_Recursive");
1227  }
1228  unsigned int b_index = fine_index%2; //This takes care that if fine_index ==2 then we select the
1229  //zeros entry in the element vector
1230 
1231  if(need_boundary_integrals && element[b_index]->at_boundary())
1232  {
1233  for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1234  {
1235  if (element[b_index]->face(face)->at_boundary()
1236  &&
1237  (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1238  element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
1239  {
1240  fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1241  pde.BoundaryRhs(fdc,local_vector,1.);
1242  }
1243  }
1244  }
1245  // if(need_faces)
1246  // {
1247  // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1248  // {
1249  // if (element[0]->neighbor_index(face) != -1)
1250  // {
1251  // fdc.ReInit(face);
1252  // pde.FaceRhs(fdc, local_vector,1.);
1253  // }
1254  // }
1255  // }
1256  if(coarse_index == 0) //Need to transfer the computed residual to the one on the coarse element
1257  {
1258  dealii::Vector<SCALAR> tmp(dofs_per_element);
1259  prolong_matrix.Tvmult(tmp,local_vector);
1260 
1261  //LocalToGlobal
1262  element[0]->get_dof_indices(local_dof_indices);
1263  for (unsigned int i = 0; i < dofs_per_element; ++i)
1264  {
1265  residual(local_dof_indices[i]) += tmp(i);
1266  }
1267  }
1268  else //Testfunctions are already the right ones...
1269  {
1270  //LocalToGlobal
1271  element[0]->get_dof_indices(local_dof_indices);
1272  for (unsigned int i = 0; i < dofs_per_element; ++i)
1273  {
1274  residual(local_dof_indices[i]) += local_vector(i);
1275  }
1276  }
1277 
1278  }//Endof the case on the finest level
1279  else
1280  {
1281  assert(fine_index != coarse_index);
1282  assert(element[fine_index]->has_children());
1283  assert(!element[coarse_index]->has_children());
1284  assert(tria_element[fine_index]->has_children());
1285  assert(!tria_element[coarse_index]->has_children());
1286 
1287  unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
1288 
1289  typename DH<dim, dim>::cell_iterator dofh_fine = element[fine_index];
1290  typename dealii::Triangulation<dim>::cell_iterator tria_fine = tria_element[fine_index];
1291 
1292  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1293  {
1294  FullMatrix<SCALAR> new_matrix(local_n_dofs);
1295  element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1296  prolong_matrix);
1297  element[fine_index] = dofh_fine->child(child);
1298  tria_element[fine_index] = tria_fine->child(child);
1299 
1300  ComputeNonlinearRhs_Recursive(pde,residual,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
1301  }
1302  }
1303  }
1304 
1305  /*******************************************************************************************/
1306 
1307  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
1308  int dim>
1309  template<typename PROBLEM, typename MATRIX, template<int, int> class DH>
1310  void
1311  IntegratorMultiMesh<INTEGRATORDATACONT, VECTOR, SCALAR, dim>::ComputeMatrix_Recursive(
1312  PROBLEM& pde, MATRIX &matrix, typename std::vector<typename DH<dim, dim>::cell_iterator>& element,
1313  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element,
1314  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
1315  Multimesh_ElementDataContainer<DH, VECTOR, dim>& edc,
1316  Multimesh_FaceDataContainer<DH, VECTOR, dim>& fdc)
1317  {
1318 
1319  if(!element[0]->has_children() && ! element[1]->has_children())
1320  {
1321  unsigned int dofs_per_element;
1322  std::vector<unsigned int> local_dof_indices;
1323 
1324  bool need_faces = pde.HasFaces();
1325  bool need_interfaces = pde.HasInterfaces();
1326  std::vector<unsigned int> boundary_equation_colors = pde.GetBoundaryEquationColors();
1327  bool need_boundary_integrals = (boundary_equation_colors.size() > 0);
1328 
1329  edc.ReInit(coarse_index,fine_index,prolong_matrix);
1330  dofs_per_element = element[0]->get_fe().dofs_per_cell;
1331 
1332  dealii::FullMatrix<SCALAR> local_matrix(dofs_per_element,
1333  dofs_per_element);
1334  local_matrix = 0;
1335 
1336  local_dof_indices.resize(0);
1337  local_dof_indices.resize(dofs_per_element, 0);
1338  pde.ElementMatrix(edc, local_matrix);
1339 
1340  //FIXME Integrate on Faces of Element[0] that contain a fine-element face.
1341  if(need_faces || need_interfaces)
1342  {
1343  throw DOpEException(" Faces on multiple meshes not implemented yet!",
1344  "IntegratorMultiMesh::ComputeMatrix_Recursive");
1345  }
1346  unsigned int b_index = fine_index%2; //This takes care that if fine_index ==2 then we select the
1347  //zeros entry in the element vector
1348  if(need_boundary_integrals && element[b_index]->at_boundary())
1349  {
1350 
1351  for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1352  {
1353  if (element[b_index]->face(face)->at_boundary()
1354  &&
1355  (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1356  element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
1357  {
1358  fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1359  pde.BoundaryMatrix(fdc, local_matrix);
1360  }
1361  }
1362  }
1363  // if(need_faces)
1364  // {
1365  // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1366  // {
1367  // if (element[0]->neighbor_index(face) != -1)
1368  // {
1369  // fdc.ReInit(face);
1370  // pde.FaceMatrix(fdc, local_matrix);
1371  // }
1372  // }
1373  // }
1374  // if( need_interfaces)
1375  // {
1376  // for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1377  // {
1378  // fdc.ReInit(face);
1379  // if (element[0]->neighbor_index(face) != -1
1380  // &&
1381  // fdc.GetMaterialId()!= fdc.GetNbrMaterialId())
1382  // {
1383  // fdc.ReInitNbr();
1384  //
1385  // nbr_dofs_per_element = fdc.GetNbrNDoFsPerElement();
1386  // nbr_local_dof_indices.resize(0);
1387  // nbr_local_dof_indices.resize(nbr_dofs_per_element, 0);
1388  // dealii::FullMatrix<SCALAR> local_interface_matrix(dofs_per_element,nbr_dofs_per_element );
1389  // local_interface_matrix = 0;
1390  //
1391  // pde.InterfaceMatrix(fdc, local_interface_matrix);
1392  //
1393  // element[0]->get_dof_indices(local_dof_indices);
1394  // element[0]->neighbor(face)->get_dof_indices(nbr_local_dof_indices);
1395  //
1396  // for (unsigned int i = 0; i < dofs_per_element; ++i)
1397  // {
1398  // for (unsigned int j = 0; j < nbr_dofs_per_element; ++j)
1399  // {
1400  // matrix.add(local_dof_indices[i], nbr_local_dof_indices[j],
1401  // local_interface_matrix(i, j));
1402  // } //endfor j
1403  // } //endfor i
1404  // }
1405  // }
1406  // }
1407  if(coarse_index == 0) //Need to transfer the computed residual to the one on the coarse element
1408  {
1409  dealii::FullMatrix<SCALAR> tmp(dofs_per_element);
1410  tmp = 0.;
1411  prolong_matrix.Tmmult(tmp,local_matrix);
1412  local_matrix = 0.;
1413  tmp.mmult(local_matrix,prolong_matrix);
1414 
1415  //LocalToGlobal
1416  element[0]->get_dof_indices(local_dof_indices);
1417  for (unsigned int i = 0; i < dofs_per_element; ++i)
1418  {
1419  for (unsigned int j = 0; j < dofs_per_element; ++j)
1420  {
1421  matrix.add(local_dof_indices[i], local_dof_indices[j],
1422  local_matrix(i, j));
1423  }
1424  }
1425  }
1426  else
1427  {
1428  //LocalToGlobal
1429  element[0]->get_dof_indices(local_dof_indices);
1430  for (unsigned int i = 0; i < dofs_per_element; ++i)
1431  {
1432  for (unsigned int j = 0; j < dofs_per_element; ++j)
1433  {
1434  matrix.add(local_dof_indices[i], local_dof_indices[j],
1435  local_matrix(i, j));
1436  }
1437  }
1438  }
1439 
1440  }//Endof the case on the finest level
1441  else
1442  {
1443  assert(fine_index != coarse_index);
1444  assert(element[fine_index]->has_children());
1445  assert(!element[coarse_index]->has_children());
1446  assert(tria_element[fine_index]->has_children());
1447  assert(!tria_element[coarse_index]->has_children());
1448 
1449  unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
1450 
1451  typename DH<dim, dim>::cell_iterator dofh_fine = element[fine_index];
1452  typename dealii::Triangulation<dim>::cell_iterator tria_fine = tria_element[fine_index];
1453 
1454  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1455  {
1456  FullMatrix<SCALAR> new_matrix(local_n_dofs);
1457  element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1458  prolong_matrix);
1459  element[fine_index] = dofh_fine->child(child);
1460  tria_element[fine_index] = tria_fine->child(child);
1461 
1462  ComputeMatrix_Recursive(pde,matrix,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
1463  }
1464  }
1465  }
1466 
1467  /*******************************************************************************************/
1468 
1469  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
1470  template<typename PROBLEM, template<int, int> class DH>
1471  SCALAR IntegratorMultiMesh<INTEGRATORDATACONT, VECTOR, SCALAR, dim>::ComputeDomainScalar_Recursive(
1472  PROBLEM& pde,
1473  typename std::vector<typename DH<dim, dim>::cell_iterator> &element,
1474  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator> &tria_element,
1475  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
1476  Multimesh_ElementDataContainer<DH, VECTOR, dim>& edc)
1477  {
1478  if(!element[0]->has_children() && ! element[1]->has_children())
1479  {
1480  SCALAR ret = 0.;
1481  edc.ReInit(coarse_index,fine_index,prolong_matrix);
1482  ret += pde.ElementFunctional(edc);
1483  return ret;
1484  } //Endof the case on the finest level
1485  else
1486  {
1487  assert(fine_index != coarse_index);
1488  assert(element[fine_index]->has_children());
1489  assert(!element[coarse_index]->has_children());
1490  assert(tria_element[fine_index]->has_children());
1491  assert(!tria_element[coarse_index]->has_children());
1492 
1493  unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
1494 
1495  typename DH<dim, dim>::cell_iterator dofh_fine = element[fine_index];
1496  typename dealii::Triangulation<dim>::cell_iterator tria_fine = tria_element[fine_index];
1497  SCALAR ret = 0.;
1498  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1499  {
1500  FullMatrix<SCALAR> new_matrix(local_n_dofs);
1501  element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1502  prolong_matrix);
1503  element[fine_index] = dofh_fine->child(child);
1504  tria_element[fine_index] = tria_fine->child(child);
1505 
1506  ret += ComputeDomainScalar_Recursive(pde,element,tria_element,new_matrix, coarse_index, fine_index, edc);
1507  }
1508  return ret;
1509  }
1510  }
1511  /*******************************************************************************************/
1512 
1513  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,int dim>
1514  template<typename PROBLEM, template<int, int> class DH>
1515  SCALAR IntegratorMultiMesh<INTEGRATORDATACONT, VECTOR, SCALAR, dim>::ComputeBoundaryScalar_Recursive(
1516  PROBLEM& pde,
1517  typename std::vector<typename DH<dim, dim>::cell_iterator> &element,
1518  typename std::vector<typename dealii::Triangulation<dim>::cell_iterator> &tria_element,
1519  const FullMatrix<SCALAR>& prolong_matrix,unsigned int coarse_index,unsigned int fine_index,
1520  Multimesh_FaceDataContainer<DH, VECTOR, dim>& fdc)
1521  {
1522  if(!element[0]->has_children() && ! element[1]->has_children())
1523  {
1524  SCALAR ret = 0.;
1525  std::vector<unsigned int> boundary_functional_colors = pde.GetBoundaryFunctionalColors();
1526  unsigned int b_index = fine_index%2; //This takes care that if fine_index ==2 then we select the
1527  //zeros entry in the element vector
1528  for (unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1529  {
1530  if (element[b_index]->face(face)->at_boundary()
1531  &&
1532  (find(boundary_functional_colors.begin(),boundary_functional_colors.end(),
1533  element[b_index]->face(face)->boundary_indicator()) != boundary_functional_colors.end()))
1534  {
1535  fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1536  ret += pde.BoundaryFunctional(fdc);
1537  }
1538  }
1539  return ret;
1540  } //Endof the case on the finest level
1541  else
1542  {
1543  assert(fine_index != coarse_index);
1544  assert(element[fine_index]->has_children());
1545  assert(!element[coarse_index]->has_children());
1546  assert(tria_element[fine_index]->has_children());
1547  assert(!tria_element[coarse_index]->has_children());
1548 
1549  unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
1550 
1551  typename DH<dim, dim>::cell_iterator dofh_fine = element[fine_index];
1552  typename dealii::Triangulation<dim>::cell_iterator tria_fine = tria_element[fine_index];
1553  SCALAR ret = 0.;
1554  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1555  {
1556  FullMatrix<SCALAR> new_matrix(local_n_dofs);
1557  element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1558  prolong_matrix);
1559  element[fine_index] = dofh_fine->child(child);
1560  tria_element[fine_index] = tria_fine->child(child);
1561 
1562  ret += ComputeBoundaryScalar_Recursive(pde,element,tria_element,new_matrix, coarse_index, fine_index, fdc);
1563  }
1564  return ret;
1565  }
1566  }
1567 //ENDOF NAMESPACE DOpE
1568 }
1569 #endif
1570 
const std::map< std::string, const VECTOR * > & GetDomainData() const
Definition: integrator_multimesh.h:964
void ComputeLocalControlConstraints(PROBLEM &pde, VECTOR &constraints)
Definition: integrator_multimesh.h:806
SCALAR ComputeBoundaryScalar(PROBLEM &pde)
Definition: integrator_multimesh.h:628
Definition: integrator_multimesh.h:62
SCALAR ComputeFaceScalar(PROBLEM &pde)
Definition: integrator_multimesh.h:710
SCALAR ComputeDomainScalar(PROBLEM &pde)
Definition: integrator_multimesh.h:529
void ComputeNonlinearLhs(PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
Definition: integrator_multimesh.h:346
Definition: multimesh_elementdatacontainer.h:55
void ComputeNonlinearRhs(PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
Definition: integrator_multimesh.h:360
void ComputeNonlinearResidual(PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
Definition: integrator_multimesh.h:256
SCALAR ComputeAlgebraicScalar(PROBLEM &pde)
Definition: integrator_multimesh.h:776
void AddParamData(std::string name, const dealii::Vector< SCALAR > *new_data)
Definition: integrator_multimesh.h:974
const DOFHANDLER< dim, dim > & GetDEALDoFHandler() const
Definition: dofhandler_wrapper.h:69
SCALAR ComputePointScalar(PROBLEM &pde)
Definition: integrator_multimesh.h:610
void ApplyNewtonBoundaryValues(PROBLEM &pde, VECTOR &u)
Definition: integrator_multimesh.h:850
void ComputeMatrix(PROBLEM &pde, MATRIX &matrix)
Definition: integrator_multimesh.h:448
~IntegratorMultiMesh()
Definition: integrator_multimesh.h:235
void DeleteParamData(std::string name)
Definition: integrator_multimesh.h:992
void ReInit()
Definition: integrator_multimesh.h:245
void DeleteDomainData(std::string name)
Definition: integrator_multimesh.h:929
void ApplyInitialBoundaryValues(PROBLEM &pde, VECTOR &u)
Definition: integrator_multimesh.h:818
IntegratorMultiMesh(INTEGRATORDATACONT &idc)
Definition: integrator_multimesh.h:225
void AddDomainData(std::string name, const VECTOR *new_data)
Definition: integrator_multimesh.h:912
void ComputeNonlinearAlgebraicResidual(PROBLEM &pde, VECTOR &residual)
Definition: integrator_multimesh.h:793
Definition: multimesh_facedatacontainer.h:53
void AddPresetRightHandSide(double s, VECTOR &residual) const
Definition: integrator_multimesh.h:947
INTEGRATORDATACONT & GetIntegratorDataContainer() const
Definition: integrator_multimesh.h:1021
Definition: dopeexception.h:35
const std::map< std::string, const dealii::Vector< SCALAR > * > & GetParamData() const
Definition: integrator_multimesh.h:1011