24 #ifndef IntegratorMultiMesh_H_
25 #define IntegratorMultiMesh_H_
27 #include <lac/vector.h>
28 #include <lac/block_sparsity_pattern.h>
29 #include <lac/block_sparse_matrix.h>
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>
60 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
72 template<
typename PROBLEM>
75 bool apply_boundary_values =
true);
76 template<
typename PROBLEM>
79 bool apply_boundary_values =
true);
80 template<
typename PROBLEM>
83 bool apply_boundary_values =
true);
84 template<
typename PROBLEM,
typename MATRIX>
88 template<
typename PROBLEM>
90 template<
typename PROBLEM>
93 template<
typename PROBLEM>
96 template<
typename PROBLEM>
99 template<
typename PROBLEM>
102 template<
typename PROBLEM>
105 template<
typename PROBLEM>
109 template<
typename PROBLEM>
112 template<
typename PROBLEM>
115 template<
typename PROBLEM,
typename MATRIX>
126 AddParamData(std::string name,
const dealii::Vector<SCALAR>* new_data);
131 inline const std::map<std::string, const dealii::Vector<SCALAR>*>&
134 inline INTEGRATORDATACONT&
136 inline const std::map<std::string, const VECTOR*>&
141 template<
template<
int,
int>
class DH>
143 InterpolateBoundaryValues(
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;
153 template<
typename PROBLEM,
template<
int,
int>
class DH>
154 inline void ComputeNonlinearResidual_Recursive(
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,
167 template<
typename PROBLEM,
template<
int,
int>
class DH>
168 inline void ComputeNonlinearRhs_Recursive(
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,
181 template<
typename PROBLEM,
typename MATRIX,
template<
int,
int>
class DH>
182 inline void ComputeMatrix_Recursive(
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,
195 template<
typename PROBLEM,
template<
int,
int>
class DH>
196 inline SCALAR ComputeDomainScalar_Recursive(
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,
207 template<
typename PROBLEM,
template<
int,
int>
class DH>
208 inline SCALAR ComputeBoundaryScalar_Recursive(
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,
215 INTEGRATORDATACONT & idc_;
217 std::map<std::string, const VECTOR*> domain_data_;
218 std::map<std::string, const dealii::Vector<SCALAR>*> param_data_;
223 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
226 INTEGRATORDATACONT& idc) :
233 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
242 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
252 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
254 template<
typename PROBLEM>
257 PROBLEM& pde, VECTOR &residual,
bool apply_boundary_values)
261 const auto& dof_handler =
262 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
264 assert(dof_handler.size() == 2);
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());
271 auto element_iter = element_list.begin();
272 auto tria_element_iter = tria_element_list.begin();
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;
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;
285 idc_.InitializeMMEDC(pde.GetUpdateFlags(),
286 *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
287 this->GetParamData(), this->GetDomainData());
288 auto& edc = idc_.GetMultimeshElementDataContainer();
290 bool need_interfaces = pde.HasInterfaces();
291 idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
292 *(pde.GetBaseProblem().GetSpaceTimeHandler()),
295 this->GetParamData(),
296 this->GetDomainData(),
298 auto& fdc = idc_.GetMultimeshFaceDataContainer();
300 for(; element_iter != element_list.end(); element_iter++)
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;
308 if(element[0]->has_children())
310 prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
316 if(element[1]->has_children())
318 prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
324 assert(element.size() ==2);
325 coarse_index = fine_index = 2;
328 ComputeNonlinearResidual_Recursive(pde,residual,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
332 AddPresetRightHandSide(-1.,residual);
334 if (apply_boundary_values)
336 ApplyNewtonBoundaryValues(pde, residual);
342 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
344 template<
typename PROBLEM>
347 PROBLEM& pde, VECTOR &residual,
bool apply_boundary_values)
350 throw DOpEException(
"This function needs to be implemented!",
"IntegratorMultiMesh::ComputeNonlinearLhs");
356 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
358 template<
typename PROBLEM>
361 PROBLEM& pde, VECTOR &residual,
bool apply_boundary_values)
365 const auto& dof_handler =
366 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
368 assert(dof_handler.size() == 2);
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());
375 auto element_iter = element_list.begin();
376 auto tria_element_iter = tria_element_list.begin();
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;
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;
389 idc_.InitializeMMEDC(pde.GetUpdateFlags(),
390 *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
391 this->GetParamData(), this->GetDomainData());
392 auto& edc = idc_.GetMultimeshElementDataContainer();
394 idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
395 *(pde.GetBaseProblem().GetSpaceTimeHandler()),
398 this->GetParamData(),
399 this->GetDomainData());
400 auto& fdc = idc_.GetMultimeshFaceDataContainer();
402 for(; element_iter != element_list.end(); element_iter++)
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;
410 if(element[0]->has_children())
412 prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
418 if(element[1]->has_children())
420 prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
426 assert(element.size() ==2);
427 coarse_index = fine_index = 2;
430 ComputeNonlinearRhs_Recursive(pde,residual,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
434 AddPresetRightHandSide(1.,residual);
436 if (apply_boundary_values)
438 ApplyNewtonBoundaryValues(pde, residual);
444 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
446 template<
typename PROBLEM,
typename MATRIX>
449 PROBLEM& pde, MATRIX &matrix)
453 const auto& dof_handler =
454 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
456 assert(dof_handler.size() == 2);
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());
463 auto element_iter = element_list.begin();
464 auto tria_element_iter = tria_element_list.begin();
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;
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;
477 idc_.InitializeMMEDC(pde.GetUpdateFlags(),
478 *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
479 this->GetParamData(), this->GetDomainData());
480 auto& edc = idc_.GetMultimeshElementDataContainer();
482 bool need_interfaces = pde.HasInterfaces();
483 idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
484 *(pde.GetBaseProblem().GetSpaceTimeHandler()),
487 this->GetParamData(),
488 this->GetDomainData(),
490 auto & fdc = idc_.GetMultimeshFaceDataContainer();
492 for(; element_iter != element_list.end(); element_iter++)
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;
500 if(element[0]->has_children())
502 prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
508 if(element[1]->has_children())
510 prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
516 assert(element.size() ==2);
517 coarse_index = fine_index = 2;
520 ComputeMatrix_Recursive(pde,matrix,element,tria_element,prolong_matrix,coarse_index,fine_index,edc,fdc);
525 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
527 template<
typename PROBLEM>
534 const auto& dof_handler =
535 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
537 assert(dof_handler.size() == 2);
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());
544 auto element_iter = element_list.begin();
545 auto tria_element_iter = tria_element_list.begin();
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;
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;
559 throw DOpEException(
"This function should not be called when faces are needed!",
560 "IntegratorMultiMesh::ComputeDomainScalar");
564 idc_.InitializeMMEDC(pde.GetUpdateFlags(),
565 *(pde.GetBaseProblem().GetSpaceTimeHandler()), element, tria_element,
566 this->GetParamData(), this->GetDomainData());
567 auto& edc = idc_.GetMultimeshElementDataContainer();
569 for(; element_iter != element_list.end(); element_iter++)
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;
577 if(element[0]->has_children())
579 prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
585 if(element[1]->has_children())
587 prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
593 assert(element.size() ==2);
594 coarse_index = fine_index = 2;
597 ret += ComputeDomainScalar_Recursive(pde,element,tria_element,prolong_matrix,coarse_index,fine_index,edc);
606 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
608 template<
typename PROBLEM>
616 ret += pde.PointFunctional(this->GetParamData(),
617 this->GetDomainData());
624 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
626 template<
typename PROBLEM>
634 const auto& dof_handler =
635 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
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());
642 auto element_iter = element_list.begin();
643 auto tria_element_iter = tria_element_list.begin();
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;
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;
656 idc_.InitializeMMFDC(pde.GetFaceUpdateFlags(),
657 *(pde.GetBaseProblem().GetSpaceTimeHandler()),
658 element, tria_element,
659 this->GetParamData(),
660 this->GetDomainData());
661 auto & fdc = idc_.GetMultimeshFaceDataContainer();
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)
667 throw DOpEException(
"No boundary colors given!",
"IntegratorMultiMesh::ComputeBoundaryScalar");
670 for(; element_iter != element_list.end(); element_iter++)
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;
678 if(element[0]->has_children())
680 prolong_matrix = IdentityMatrix(element[1]->get_fe().dofs_per_cell);
686 if(element[1]->has_children())
688 prolong_matrix = IdentityMatrix(element[0]->get_fe().dofs_per_cell);
694 assert(element.size() ==2);
695 coarse_index = fine_index = 2;
698 ret += ComputeBoundaryScalar_Recursive(pde,element,tria_element,prolong_matrix,coarse_index,fine_index,fdc);
706 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
708 template<
typename PROBLEM>
714 throw DOpEException(
"This function needs to be implemented!",
"IntegratorMultiMesh::ComputeFaceScalar");
772 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
774 template<
typename PROBLEM>
782 ret = pde.AlgebraicFunctional(this->GetParamData(),
783 this->GetDomainData());
790 template <
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
int dim>
791 template<
typename PROBLEM>
796 pde.AlgebraicResidual(residual,this->GetParamData(),this->GetDomainData());
798 AddPresetRightHandSide(-1.,residual);
803 template <
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
int dim>
804 template<
typename PROBLEM>
809 pde.ComputeLocalControlConstraints(constraints,this->GetParamData(),this->GetDomainData());
814 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
816 template<
typename PROBLEM>
819 PROBLEM& pde, VECTOR &u)
826 std::vector<unsigned int> dirichlet_colors = pde.GetDirichletColors();
827 for (
unsigned int i = 0; i < dirichlet_colors.size(); i++)
829 unsigned int color = dirichlet_colors[i];
830 std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
831 std::map<unsigned int, SCALAR> boundary_values;
833 InterpolateBoundaryValues(pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0], color, pde.GetDirichletValues(color, this->GetParamData(),
834 this->GetDomainData()), boundary_values, comp_mask );
836 for (
typename std::map<unsigned int, SCALAR>::const_iterator p =
837 boundary_values.begin(); p != boundary_values.end(); p++)
839 u(p->first) = p->second;
846 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
848 template<
typename PROBLEM>
851 PROBLEM& pde, VECTOR &u)
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++)
860 unsigned int color = dirichlet_colors[i];
861 std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
862 std::map<unsigned int, SCALAR> boundary_values;
865 InterpolateBoundaryValues(
866 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0],
867 color, dealii::ZeroFunction<dim>(comp_mask.size()),
868 boundary_values, comp_mask);
870 for (
typename std::map<unsigned int, SCALAR>::const_iterator p =
871 boundary_values.begin(); p != boundary_values.end(); p++)
873 u(p->first) = p->second;
879 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
881 template<
typename PROBLEM,
typename MATRIX>
884 PROBLEM& pde, MATRIX& matrix, VECTOR &rhs, VECTOR &sol)
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++)
893 unsigned int color = dirichlet_colors[i];
894 std::vector<bool> comp_mask = pde.GetDirichletCompMask(color);
895 std::map<unsigned int, SCALAR> boundary_values;
897 InterpolateBoundaryValues(
898 pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0],
899 color, dealii::ZeroFunction<dim>(comp_mask.size()),
900 boundary_values, comp_mask);
902 dealii::MatrixTools::apply_boundary_values(boundary_values, matrix,
909 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
913 std::string name,
const VECTOR* new_data)
915 if (domain_data_.find(name) != domain_data_.end())
918 "Adding multiple Data with name " + name +
" is prohibited!",
919 "IntegratorMultiMesh::AddDomainData");
921 domain_data_.insert(std::pair<std::string, const VECTOR*>(name, new_data));
926 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
932 typename std::map<std::string, const VECTOR *>::iterator it =
933 domain_data_.find(name);
934 if (it == domain_data_.end())
937 "Deleting Data " + name +
" is impossible! Data not found",
938 "IntegratorMultiMesh::DeleteDomainData");
940 domain_data_.erase(it);
944 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
948 VECTOR &residual)
const
950 typename std::map<std::string, const VECTOR *>::const_iterator it =
951 domain_data_.find(
"fixed_rhs");
952 if (it != domain_data_.end())
954 assert(residual.size() == it->second->size());
955 residual.add(s,*(it->second));
961 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
963 const std::map<std::string, const VECTOR*>&
971 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
975 std::string name,
const dealii::Vector<SCALAR>* new_data)
977 if (param_data_.find(name) != param_data_.end())
980 "Adding multiple Data with name " + name +
" is prohibited!",
981 "IntegratorMultiMesh::AddParamData");
984 std::pair<std::string,
const dealii::Vector<SCALAR>*>(name, new_data));
989 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
995 typename std::map<std::string, const dealii::Vector<SCALAR>*>::iterator
996 it = param_data_.find(name);
997 if (it == param_data_.end())
1000 "Deleting Data " + name +
" is impossible! Data not found",
1001 "IntegratorMultiMesh::DeleteParamData");
1003 param_data_.erase(it);
1008 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
1010 const std::map<std::string, const dealii::Vector<SCALAR>*>&
1018 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
1029 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
1031 template<
template<
int,
int>
class DH>
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
1039 dealii::VectorTools::interpolate_boundary_values(
1041 boundary_values, comp_mask);
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(
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)
1056 if(!element[0]->has_children() && ! element[1]->has_children())
1058 unsigned int dofs_per_element;
1059 dealii::Vector<SCALAR> local_vector;
1060 std::vector<unsigned int> local_dof_indices;
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();
1067 edc.
ReInit(coarse_index,fine_index,prolong_matrix);
1069 dofs_per_element = element[0]->get_fe().dofs_per_cell;
1071 local_vector.reinit(dofs_per_element);
1074 local_dof_indices.resize(0);
1075 local_dof_indices.resize(dofs_per_element, 0);
1079 pde.ElementEquation(edc, local_vector, 1., 1.);
1080 pde.ElementRhs(edc, local_vector, -1.);
1083 if(need_faces || need_interfaces)
1085 throw DOpEException(
" Faces on multiple meshes not implemented yet!",
1086 "IntegratorMultiMesh::ComputeNonlinearResidual_Recursive");
1088 unsigned int b_index = fine_index%2;
1091 if(need_boundary_integrals && element[b_index]->at_boundary())
1093 for (
unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1095 if (element[b_index]->face(face)->at_boundary()
1097 (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1098 element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
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.);
1132 if(coarse_index == 0)
1134 dealii::Vector<SCALAR> tmp(dofs_per_element);
1135 prolong_matrix.Tvmult(tmp,local_vector);
1138 element[0]->get_dof_indices(local_dof_indices);
1139 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1141 residual(local_dof_indices[i]) += tmp(i);
1147 element[0]->get_dof_indices(local_dof_indices);
1148 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1150 residual(local_dof_indices[i]) += local_vector(i);
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());
1163 unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
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];
1168 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1170 FullMatrix<SCALAR> new_matrix(local_n_dofs);
1171 element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1173 element[fine_index] = dofh_fine->child(child);
1174 tria_element[fine_index] = tria_fine->child(child);
1176 ComputeNonlinearResidual_Recursive(pde,residual,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
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(
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)
1192 if(!element[0]->has_children() && ! element[1]->has_children())
1194 unsigned int dofs_per_element;
1195 dealii::Vector<SCALAR> local_vector;
1196 std::vector<unsigned int> local_dof_indices;
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();
1203 edc.ReInit(coarse_index,fine_index,prolong_matrix);
1205 dofs_per_element = element[0]->get_fe().dofs_per_cell;
1207 local_vector.reinit(dofs_per_element);
1210 local_dof_indices.resize(0);
1211 local_dof_indices.resize(dofs_per_element, 0);
1215 pde.ElementRhs(edc, local_vector, 1.);
1220 throw DOpEException(
" Faces on multiple meshes not implemented yet!",
1221 "IntegratorMultiMesh::ComputeNonlinearRhs_Recursive");
1223 if(need_interfaces )
1225 throw DOpEException(
" Faces on multiple meshes not implemented yet!",
1226 "IntegratorMultiMesh::ComputeNonlinearRhs_Recursive");
1228 unsigned int b_index = fine_index%2;
1231 if(need_boundary_integrals && element[b_index]->at_boundary())
1233 for (
unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1235 if (element[b_index]->face(face)->at_boundary()
1237 (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1238 element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
1240 fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1241 pde.BoundaryRhs(fdc,local_vector,1.);
1256 if(coarse_index == 0)
1258 dealii::Vector<SCALAR> tmp(dofs_per_element);
1259 prolong_matrix.Tvmult(tmp,local_vector);
1262 element[0]->get_dof_indices(local_dof_indices);
1263 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1265 residual(local_dof_indices[i]) += tmp(i);
1271 element[0]->get_dof_indices(local_dof_indices);
1272 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1274 residual(local_dof_indices[i]) += local_vector(i);
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());
1287 unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
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];
1292 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1294 FullMatrix<SCALAR> new_matrix(local_n_dofs);
1295 element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1297 element[fine_index] = dofh_fine->child(child);
1298 tria_element[fine_index] = tria_fine->child(child);
1300 ComputeNonlinearRhs_Recursive(pde,residual,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
1307 template<
typename INTEGRATORDATACONT,
typename VECTOR,
typename SCALAR,
1309 template<
typename PROBLEM,
typename MATRIX,
template<
int,
int>
class DH>
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)
1319 if(!element[0]->has_children() && ! element[1]->has_children())
1321 unsigned int dofs_per_element;
1322 std::vector<unsigned int> local_dof_indices;
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);
1329 edc.ReInit(coarse_index,fine_index,prolong_matrix);
1330 dofs_per_element = element[0]->get_fe().dofs_per_cell;
1332 dealii::FullMatrix<SCALAR> local_matrix(dofs_per_element,
1336 local_dof_indices.resize(0);
1337 local_dof_indices.resize(dofs_per_element, 0);
1338 pde.ElementMatrix(edc, local_matrix);
1341 if(need_faces || need_interfaces)
1343 throw DOpEException(
" Faces on multiple meshes not implemented yet!",
1344 "IntegratorMultiMesh::ComputeMatrix_Recursive");
1346 unsigned int b_index = fine_index%2;
1348 if(need_boundary_integrals && element[b_index]->at_boundary())
1351 for (
unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1353 if (element[b_index]->face(face)->at_boundary()
1355 (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),
1356 element[b_index]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
1358 fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1359 pde.BoundaryMatrix(fdc, local_matrix);
1407 if(coarse_index == 0)
1409 dealii::FullMatrix<SCALAR> tmp(dofs_per_element);
1411 prolong_matrix.Tmmult(tmp,local_matrix);
1413 tmp.mmult(local_matrix,prolong_matrix);
1416 element[0]->get_dof_indices(local_dof_indices);
1417 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1419 for (
unsigned int j = 0; j < dofs_per_element; ++j)
1421 matrix.add(local_dof_indices[i], local_dof_indices[j],
1422 local_matrix(i, j));
1429 element[0]->get_dof_indices(local_dof_indices);
1430 for (
unsigned int i = 0; i < dofs_per_element; ++i)
1432 for (
unsigned int j = 0; j < dofs_per_element; ++j)
1434 matrix.add(local_dof_indices[i], local_dof_indices[j],
1435 local_matrix(i, j));
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());
1449 unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
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];
1454 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1456 FullMatrix<SCALAR> new_matrix(local_n_dofs);
1457 element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1459 element[fine_index] = dofh_fine->child(child);
1460 tria_element[fine_index] = tria_fine->child(child);
1462 ComputeMatrix_Recursive(pde,matrix,element,tria_element,new_matrix, coarse_index, fine_index, edc, fdc);
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(
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)
1478 if(!element[0]->has_children() && ! element[1]->has_children())
1481 edc.ReInit(coarse_index,fine_index,prolong_matrix);
1482 ret += pde.ElementFunctional(edc);
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());
1493 unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
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];
1498 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1500 FullMatrix<SCALAR> new_matrix(local_n_dofs);
1501 element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1503 element[fine_index] = dofh_fine->child(child);
1504 tria_element[fine_index] = tria_fine->child(child);
1506 ret += ComputeDomainScalar_Recursive(pde,element,tria_element,new_matrix, coarse_index, fine_index, edc);
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(
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)
1522 if(!element[0]->has_children() && ! element[1]->has_children())
1525 std::vector<unsigned int> boundary_functional_colors = pde.GetBoundaryFunctionalColors();
1526 unsigned int b_index = fine_index%2;
1528 for (
unsigned int face=0; face < dealii::GeometryInfo<dim>::faces_per_cell; ++face)
1530 if (element[b_index]->face(face)->at_boundary()
1532 (find(boundary_functional_colors.begin(),boundary_functional_colors.end(),
1533 element[b_index]->face(face)->boundary_indicator()) != boundary_functional_colors.end()))
1535 fdc.ReInit(coarse_index,fine_index,prolong_matrix,face);
1536 ret += pde.BoundaryFunctional(fdc);
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());
1549 unsigned int local_n_dofs = element[coarse_index]->get_fe().dofs_per_cell;
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];
1554 for (
unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell;++child)
1556 FullMatrix<SCALAR> new_matrix(local_n_dofs);
1557 element[coarse_index]->get_fe().get_prolongation_matrix(child).mmult (new_matrix,
1559 element[fine_index] = dofh_fine->child(child);
1560 tria_element[fine_index] = tria_fine->child(child);
1562 ret += ComputeBoundaryScalar_Recursive(pde,element,tria_element,new_matrix, coarse_index, fine_index, fdc);
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