DOpE
integratormixeddims.h
Go to the documentation of this file.
1 
24 #ifndef IntegratorMixed_H_
25 #define IntegratorMixed_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 
34 #include <base/function.h>
35 
36 #include <dofs/dof_tools.h>
37 
38 #include <fe/mapping_q1.h>
39 
40 #include <vector>
41 
42 #include "elementdatacontainer.h"
43 #include "facedatacontainer.h"
44 #include "optproblemcontainer.h"
45 #if DEAL_II_VERSION_GTE(7,3)
46  #include <base/types.h>
47  #endif
48 
49 namespace DOpE
50 {
68 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
69  int dimhigh>
71 {
72  public:
73  IntegratorMixedDimensions(INTEGRATORDATACONT& idc);
74 
76 
77  void ReInit();
78 
79  template<typename PROBLEM>
80  void ComputeNonlinearResidual(PROBLEM& pde, VECTOR &residual, bool apply_boundary_values = true);
81  template<typename PROBLEM, typename MATRIX>
82  void ComputeMatrix(PROBLEM& pde, MATRIX &matrix);
83  template<typename PROBLEM>
84  void
85  ComputeNonlinearRhs(PROBLEM& pde, VECTOR &residual, bool apply_boundary_values = true);
86 
87  template<typename PROBLEM>
88  void ComputeLocalControlConstraints (PROBLEM& pde, VECTOR &constraints);
89  template<typename PROBLEM>
90  SCALAR ComputeDomainScalar(PROBLEM& pde);
91  template<typename PROBLEM>
92  SCALAR ComputePointScalar(PROBLEM& pde);
93  template<typename PROBLEM>
94  SCALAR ComputeBoundaryScalar(PROBLEM& pde);
95  template<typename PROBLEM>
96  SCALAR ComputeFaceScalar(PROBLEM& pde);
97  template<typename PROBLEM>
98  SCALAR ComputeAlgebraicScalar(PROBLEM& pde);
99 
100  template<typename PROBLEM>
101  void ApplyInitialBoundaryValues(PROBLEM& pde, VECTOR &u);
102  template<typename PROBLEM>
103  void ApplyTransposedInitialBoundaryValues(PROBLEM& pde, VECTOR &u, SCALAR scale);
104  template<typename PROBLEM>
105  void ApplyNewtonBoundaryValues(PROBLEM& pde, VECTOR &u);
106  template<typename PROBLEM, typename MATRIX>
107  void ApplyNewtonBoundaryValues(PROBLEM& pde, MATRIX &matrix, VECTOR &rhs, VECTOR &sol);
108 
109  inline void AddDomainData(std::string name, const VECTOR* new_data);
110  inline void DeleteDomainData(std::string name);
111 
112  inline void AddParamData(std::string name, const dealii::Vector<SCALAR>* new_data);
113  inline void DeleteParamData(std::string name);
114 
115  inline const INTEGRATORDATACONT& GetIntegratorDataContainer() const;
116 
117  protected:
118  inline const std::map<std::string, const VECTOR*>& GetDomainData() const;
119  inline const std::map<std::string, const dealii::Vector<SCALAR>*>& GetParamData() const;
120 
121  inline void AddPresetRightHandSide(double s, dealii::Vector<SCALAR> &residual) const;
122 
123  private:
124  INTEGRATORDATACONT & idc_;
125 
126  std::map<std::string, const VECTOR*> domain_data_;
127  std::map<std::string, const dealii::Vector<SCALAR>*> param_data_;
128 };
129 
130 /**********************************Implementation*******************************************/
131 
132  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
133  int dimlow, int dimhigh>
134  IntegratorMixedDimensions<INTEGRATORDATACONT, VECTOR, SCALAR, dimlow,
135  dimhigh>::IntegratorMixedDimensions(INTEGRATORDATACONT& idc) :
136  idc_(idc)
137  {
138  }
139 
140 /**********************************Implementation*******************************************/
141 
142 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
143  int dimhigh>
145 {
146 
147 }
148 
149 /*******************************************************************************************/
150 
151 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
152  int dimhigh>
154 {
155 
156 }
157 
158 /*******************************************************************************************/
159 
160 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
161  int dimhigh>
162 template<typename PROBLEM>
164  VECTOR &residual,
165  bool apply_boundary_values)
166 {
167  {
168  residual = 0.;
169  // Begin integration
170  const unsigned int dofs =
171  pde.GetBaseProblem().GetSpaceTimeHandler()->GetControlNDoFs();
172 
173  dealii::Vector<SCALAR> local_vector(dofs);
174  const auto& dof_handler =
175  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
176  auto element = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerBeginActive();
177  auto endc = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerEnd();
178 
179  // Initialize the data containers.
180  idc_.InitializeEDC(pde.GetUpdateFlags(),
181  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element,
182  this->GetParamData(), this->GetDomainData());
183  auto& edc = idc_.GetElementDataContainer();
184 
185  bool need_faces = pde.HasFaces();
186  std::vector<unsigned int> boundary_equation_colors = pde.GetBoundaryEquationColors();
187  bool need_boundary_integrals = (boundary_equation_colors.size() > 0);
188  idc_.InitializeFDC(pde.GetFaceUpdateFlags(),
189  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
190  element,
191  this->GetParamData(),
192  this->GetDomainData());
193  auto & fdc = idc_.GetFaceDataContainer();
194 
195  for (; element[0] != endc[0]; element[0]++)
196  {
197  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
198  {
199  if (element[dh] == endc[dh])
200  {
201  throw DOpEException("Elementnumbers in DoFHandlers are not matching!",
202  "IntegratorMixedDimensions::ComputeNonlinearResidual");
203  }
204  }
205 
206  local_vector = 0;
207  edc.ReInit();
208 
209  pde.ElementRhs(edc,local_vector, -1.);
210 
211  if(need_boundary_integrals)
212  {
213  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
214  {
215  if (element[0]->face(face)->at_boundary()
216  &&
217  (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),element[0]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
218  {
219  fdc.ReInit(face);
220  pde.BoundaryRhs(fdc,local_vector,-1.);
221  }
222  }
223  }
224  if(need_faces)
225  {
226  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
227  {
228  if (element[0]->neighbor_index(face) != -1)
229  {
230  fdc.ReInit(face);
231  pde.FaceRhs(fdc,local_vector,-1.);
232  }
233  }
234  }
235  //LocalToGlobal
236  for (unsigned int i = 0; i < dofs; ++i)
237  {
238  residual(i) += local_vector(i);
239  }
240 
241  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
242  {
243  element[dh]++;
244  }
245  }
246  //The Equation should not be space dependend.
247  local_vector = 0;
248  pde.ElementEquation(edc, local_vector, 1., 1.);
249 
250  for (unsigned int i = 0; i < dofs; ++i)
251  {
252  residual(i) += local_vector(i);
253  }
254 
255  if (pde.HasControlInDirichletData())
256  {
257  ApplyTransposedInitialBoundaryValues(pde,residual, -1.);
258  }
259  //Check if some preset righthandside exists.
260  local_vector = 0;
261  AddPresetRightHandSide(-1.,local_vector);
262  for (unsigned int i = 0; i < dofs; ++i)
263  {
264  residual(i) += local_vector(i);
265  }
266 
267  if (apply_boundary_values)
268  {
269  ApplyNewtonBoundaryValues(pde,residual);
270  }
271  }
272 }
273 /*******************************************************************************************/
274 
275 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
276  int dimhigh>
277 template<typename PROBLEM>
279  VECTOR &residual,
280  bool apply_boundary_values)
281 {
282  {
283  residual = 0.;
284  // Begin integration
285  const unsigned int dofs =
286  pde.GetBaseProblem().GetSpaceTimeHandler()->GetControlNDoFs();
287 
288  dealii::Vector<SCALAR> local_vector(dofs);
289  const auto& dof_handler =
290  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
291  auto element = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerBeginActive();
292  auto endc = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerEnd();
293 
294  // Initialize the data containers.
295  idc_.InitializeEDC(pde.GetUpdateFlags(),
296  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element,
297  this->GetParamData(), this->GetDomainData());
298  auto& edc = idc_.GetElementDataContainer();
299 
300  bool need_faces = pde.HasFaces();
301  std::vector<unsigned int> boundary_equation_colors = pde.GetBoundaryEquationColors();
302  bool need_boundary_integrals = (boundary_equation_colors.size() > 0);
303  idc_.InitializeFDC(pde.GetFaceUpdateFlags(),
304  *(pde.GetBaseProblem().GetSpaceTimeHandler()),
305  element,
306  this->GetParamData(),
307  this->GetDomainData());
308  auto & fdc = idc_.GetFaceDataContainer();
309 
310  for (; element[0] != endc[0]; element[0]++)
311  {
312  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
313  {
314  if (element[dh] == endc[dh])
315  {
316  throw DOpEException("Elementnumbers in DoFHandlers are not matching!",
317  "IntegratorMixedDimensions::ComputeNonlinearResidual");
318  }
319  }
320 
321  local_vector = 0;
322  edc.ReInit();
323 
324  pde.ElementRhs(edc,local_vector, 1.);
325 
326  if(need_boundary_integrals)
327  {
328  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
329  {
330  if (element[0]->face(face)->at_boundary()
331  &&
332  (find(boundary_equation_colors.begin(),boundary_equation_colors.end(),element[0]->face(face)->boundary_indicator()) != boundary_equation_colors.end()))
333  {
334  fdc.ReInit(face);
335  pde.BoundaryRhs(fdc,local_vector,1.);
336  }
337  }
338  }
339  if(need_faces)
340  {
341  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
342  {
343  if (element[0]->neighbor_index(face) != -1)
344  {
345  fdc.ReInit(face);
346  pde.FaceRhs(fdc,local_vector,1.);
347  }
348  }
349  }
350  //LocalToGlobal
351  for (unsigned int i = 0; i < dofs; ++i)
352  {
353  residual(i) += local_vector(i);
354  }
355 
356  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
357  {
358  element[dh]++;
359  }
360  }
361 
362  if (pde.HasControlInDirichletData())
363  {
364  ApplyTransposedInitialBoundaryValues(pde,residual, -1.);
365  }
366  //Check if some preset righthandside exists.
367  local_vector = 0;
368  AddPresetRightHandSide(-1.,local_vector);
369  for (unsigned int i = 0; i < dofs; ++i)
370  {
371  residual(i) += local_vector(i);
372  }
373 
374  if (apply_boundary_values)
375  {
376  ApplyNewtonBoundaryValues(pde,residual);
377  }
378  }
379 }
380 
381 /*******************************************************************************************/
382 
383 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
384  int dimhigh>
385 template<typename PROBLEM, typename MATRIX>
387 {
388  throw DOpEException("You should not use this function, try VoidLinearSolver instead.",
389  "IntegratorMixedDimensions::ComputeMatrix");
390 }
391 
392 /*******************************************************************************************/
393 
394 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
395  int dimhigh>
396 template<typename PROBLEM>
398 {
399  constraints = 0.;
400  pde.ComputeLocalControlConstraints(constraints,this->GetParamData(),this->GetDomainData());
401 }
402 /*******************************************************************************************/
403 
404  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
405  int dimlow, int dimhigh>
406  template<typename PROBLEM>
407  SCALAR
408  IntegratorMixedDimensions<INTEGRATORDATACONT, VECTOR, SCALAR, dimlow,
409  dimhigh>::ComputeDomainScalar(PROBLEM& pde)
410  {
411  {
412  SCALAR ret = 0.;
413  const unsigned int n_q_points = this->GetQuadratureFormula()->size();
414 
415  const auto& dof_handler =
416  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
417  auto element = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerBeginActive();
418  auto endc = pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandlerEnd();
419 
420  idc_.InitializeEDC(pde.GetUpdateFlags(),
421  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element,
422  this->GetParamData(), this->GetDomainData());
423  auto& edc = idc_.GetElementDataContainer();
424 
425  for (; element[0] != endc[0]; element[0]++)
426  {
427  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
428  {
429  if (element[dh] == endc[dh])
430  {
431  throw DOpEException("Elementnumbers in DoFHandlers are not matching!",
432  "IntegratorMixedDimensions::ComputeDomainScalar");
433  }
434  }
435 
436  ret += pde.ElementFunctional(edc);
437 
438  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
439  {
440  element[dh]++;
441  }
442  }
443  return ret;
444  }
445 }
446 /*******************************************************************************************/
447 
448 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
449  int dimhigh>
450 template<typename PROBLEM>
452 {
453  if (pde.GetFEValuesNeededToBeInitialized())
454  {
455  this->InitializeFEValues();
456  }
457 
458  {
459  SCALAR ret = 0.;
460 
461  ret += pde.PointFunctional(this->GetParamData(), this->GetDomainData());
462 
463  return ret;
464  }
465 }
466 /*******************************************************************************************/
467 
468 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
469  int dimhigh>
470 template<typename PROBLEM>
472 {
473 
474  {
475  SCALAR ret = 0.;
476  // Begin integration
477 
478  const std::vector<const DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler >*>& dof_handler =
479  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
480  std::vector<typename DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler >::active_element_iterator>
481  element(dof_handler.size());
482  std::vector<typename DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler>::active_element_iterator>
483  endc(dof_handler.size());
484 
485  for (unsigned int dh = 0; dh < dof_handler.size(); dh++)
486  {
487  element[dh] = dof_handler[dh]->begin_active();
488  endc[dh] = dof_handler[dh]->end();
489  }
490 
491  // Generate the data containers.
492  FaceDataContainer<dealii::DoFHandler, VECTOR, dimhigh> fdc(*(this->GetFaceQuadratureFormula()),
493  pde.GetFaceUpdateFlags(),
494  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element,
495  this->GetParamData(), this->GetDomainData());
496 
497  std::vector<unsigned int> boundary_functional_colors =
498  pde.GetBoundaryFunctionalColors();
499  bool need_boundary_integrals = (boundary_functional_colors.size() > 0);
500  if (!need_boundary_integrals)
501  {
502  throw DOpEException("No boundary colors given!",
503  "IntegratorMixedDimensions::ComputeBoundaryScalar");
504  }
505 
506  for (; element[0] != endc[0]; element[0]++)
507  {
508  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
509  {
510  if (element[dh] == endc[dh])
511  {
512  throw DOpEException("Elementnumbers in DoFHandlers are not matching!",
513  "IntegratorMixedDimensions::ComputeBoundaryScalar");
514  }
515  }
516 
517  if(need_boundary_integrals)
518  {
519  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
520  {
521  if (element[0]->face(face)->at_boundary()
522  &&
523  (find(boundary_functional_colors.begin(),boundary_functional_colors.end(),
524  element[0]->face(face)->boundary_indicator()) != boundary_functional_colors.end()))
525  {
526 // pde.GetBaseProblem().GetSpaceTimeHandler()->ComputeFaceFEValues(element, face, pde.GetType());
527  fdc.ReInit(face);
528  ret += pde.BoundaryFunctional(fdc);
529  }
530  }
531  }
532  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
533  {
534  element[dh]++;
535  }
536  }
537  return ret;
538  }
539 }
540 /*******************************************************************************************/
541 
542 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
543  int dimhigh>
544 template<typename PROBLEM>
546 {
547  {
548  SCALAR ret = 0.;
549 
550  // Begin integration
551  const std::vector<const DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler >*>& dof_handler =
552  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler();
553  std::vector<typename DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler>::active_element_iterator>
554  element(dof_handler.size());
555  std::vector<typename DOpEWrapper::DoFHandler<dimhigh, dealii::DoFHandler>::active_element_iterator>
556  endc(dof_handler.size());
557 
558  for (unsigned int dh = 0; dh < dof_handler.size(); dh++)
559  {
560  element[dh] = dof_handler[dh]->begin_active();
561  endc[dh] = dof_handler[dh]->end();
562  }
563  // Generate the data containers.
564  FaceDataContainer<dealii::DoFHandler, VECTOR, dimhigh> fdc(*(this->GetFaceQuadratureFormula()),
565  pde.GetFaceUpdateFlags(),
566  *(pde.GetBaseProblem().GetSpaceTimeHandler()), element,
567  this->GetParamData(), this->GetDomainData());
568 
569  bool need_faces = pde.HasFaces();
570  if (!need_faces)
571  {
572  throw DOpEException("No faces required!", "IntegratorMixedDimensions::ComputeFaceScalar");
573  }
574 
575  for (; element[0] != endc[0]; element[0]++)
576  {
577  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
578  {
579  if (element[dh] == endc[dh])
580  {
581  throw DOpEException("Elementnumbers in DoFHandlers are not matching!",
582  "IntegratorMixedDimensions::ComputeFaceScalar");
583  }
584  }
585 
586  if(need_faces)
587  {
588  for (unsigned int face=0; face < dealii::GeometryInfo<dimhigh>::faces_per_cell; ++face)
589  {
590  fdc.ReInit(face);
591  ret +=pde.FaceFunctional(fdc);
592  }
593  }
594  for (unsigned int dh = 1; dh < dof_handler.size(); dh++)
595  {
596  element[dh]++;
597  }
598  }
599  return ret;
600  }
601 }
602 
603 /*******************************************************************************************/
604 
605 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
606  int dimhigh>
607 template<typename PROBLEM>
609 {
610  SCALAR ret = 0.;
611  ret = pde.AlgebraicFunctional(this->GetParamData(), this->GetDomainData());
612  return ret;
613 }
614 /*******************************************************************************************/
615 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
616  int dimhigh>
617 template<typename PROBLEM>
619  VECTOR &/*u*/)
620 {
621 }
622 /*******************************************************************************************/
623 
624 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
625  int dimhigh>
626 template<typename PROBLEM>
629  VECTOR &u,
630  SCALAR scale)
631 {
632  //Das macht nur sinn, wenn es um "Transponierte Dirichletdaten geht.
633  unsigned int dofs =
634  pde.GetBaseProblem().GetSpaceTimeHandler()->GetControlNDoFs();
635  dealii::Vector<SCALAR> local_vector(dofs);
636 
637  std::vector<unsigned int> dirichlet_colors = pde.GetTransposedDirichletColors();
638  std::vector<bool> selected_components;
639  if (dirichlet_colors.size() > 0)
640  {
641  selected_components.resize(
642  pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0]->n_dofs());
643  const std::vector<Point<dimhigh> >& support_points =
644  pde.GetBaseProblem().GetSpaceTimeHandler()->GetMapDoFToSupportPoints();
645 
646  for (unsigned int i = 0; i < dirichlet_colors.size(); i++)
647  {
648  unsigned int color = dirichlet_colors[i];
649  std::vector<bool> comp_mask = pde.GetTransposedDirichletCompMask(color);
650  std::vector<bool> current_comp(comp_mask.size(), false);
651 #if DEAL_II_VERSION_GTE(7,3)
652  std::set<types::boundary_id> boundary_indicators;
653 #else
654  std::set<unsigned char> boundary_indicators;
655 #endif
656  boundary_indicators.insert(color);
657  for (unsigned int j = 0; j < comp_mask.size(); j++)
658  {
659  if (j > 0)
660  current_comp[j - 1] = false;
661  if (comp_mask[j])
662  {
663  current_comp[j] = true;
664  //Hole eine Liste der DoFs auf dem Rand und die zugehoerigen Knoten
665  DoFTools::extract_boundary_dofs(pde.GetBaseProblem().GetSpaceTimeHandler()->GetDoFHandler()[0]->GetDEALDoFHandler(),
666  current_comp, selected_components, boundary_indicators);
667  }
669  pde.GetTransposedDirichletValues(color,
670  this->GetParamData(),
671  this->GetDomainData());
672  for (unsigned int k = 0; k < selected_components.size(); k++)
673  {
674  if (selected_components[k])
675  {
676  local_vector = 0.;
677  DD.value(support_points[k], j, k, local_vector);
678  for (unsigned int l = 0; l < dofs; ++l)
679  {
680  u(l) += scale * local_vector(l);
681  }
682  }
683  }
684  }
685  //end loop over components
686  }
687  }
688 }
689 
690 /*******************************************************************************************/
691 
692 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
693  int dimhigh>
694 template<typename PROBLEM>
696  VECTOR &/*u*/)
697 {
698  //We don't need those in the mixed case...
699 }
700 /*******************************************************************************************/
701 
702 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
703  int dimhigh>
704 template<typename PROBLEM, typename MATRIX>
706  MATRIX &/*matrix*/,
707  VECTOR &/*rhs*/,
708  VECTOR &/*sol*/)
709 {
710  //We don't need those in the mixed case...
711 }
712 
713 /*******************************************************************************************/
714 
715 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
716  int dimhigh>
718  std::string name,
719  const VECTOR* new_data)
720 {
721  if (domain_data_.find(name) != domain_data_.end())
722  {
723  throw DOpEException("Adding multiple Data with name " + name + " is prohibited!",
724  "IntegratorMixedDimensions::AddDomainData");
725  }
726  domain_data_.insert(std::pair<std::string, const VECTOR*>(name, new_data));
727 }
728 
729 /*******************************************************************************************/
730 
731 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
732  int dimhigh>
734  std::string name)
735 {
736  typename std::map<std::string, const VECTOR *>::iterator it = domain_data_.find(name);
737  if (it == domain_data_.end())
738  {
739  throw DOpEException("Deleting Data " + name + " is impossible! Data not found",
740  "IntegratorMixedDimensions::DeleteDomainData");
741  }
742  domain_data_.erase(it);
743 }
744 
745 /*******************************************************************************************/
746 
747 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
748  int dimhigh>
749 const std::map<std::string, const VECTOR*>& IntegratorMixedDimensions<INTEGRATORDATACONT, VECTOR,
750  SCALAR, dimlow, dimhigh>::GetDomainData() const
751 {
752  return domain_data_;
753 }
754 
755 /*******************************************************************************************/
756 
757 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
758  int dimhigh>
760  std::string name,
761  const dealii::Vector<
762  SCALAR>* new_data)
763 {
764  if (param_data_.find(name) != param_data_.end())
765  {
766  throw DOpEException("Adding multiple Data with name " + name + " is prohibited!",
767  "IntegratorMixedDimensions::AddParamData");
768  }
769  param_data_.insert(std::pair<std::string, const dealii::Vector<SCALAR>*>(name, new_data));
770 }
771 
772 /*******************************************************************************************/
773 
774 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
775  int dimhigh>
777  std::string name)
778 {
779  typename std::map<std::string, const dealii::Vector<SCALAR>*>::iterator it =
780  param_data_.find(name);
781  if (it == param_data_.end())
782  {
783  throw DOpEException("Deleting Data " + name + " is impossible! Data not found",
784  "IntegratorMixedDimensions::DeleteParamData");
785  }
786  param_data_.erase(it);
787 }
788 
789 /*******************************************************************************************/
790 
791 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dimlow,
792  int dimhigh>
793 const std::map<std::string, const dealii::Vector<SCALAR>*>& IntegratorMixedDimensions<
794  INTEGRATORDATACONT, VECTOR, SCALAR, dimlow, dimhigh>::GetParamData() const
795 {
796  return param_data_;
797 }
798 /*******************************************************************************************/
799 
800 template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
801  int dimlow, int dimhigh>
802  const INTEGRATORDATACONT&
804  {
805  return idc_;
806  }
807 
808  /*******************************************************************************************/
809 
810  template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR,
811  int dimlow, int dimhigh>
812  void
814  dealii::Vector<SCALAR> &residual) const
815  {
816  typename std::map<std::string, const dealii::Vector<SCALAR>*>::const_iterator it =
817  param_data_.find("fixed_rhs");
818  if (it != param_data_.end())
819  {
820  assert(residual.size() == it->second->size());
821  residual.add(s,*(it->second));
822  }
823  }
824 
825 /*******************************************************************************************/
826 }
827 #endif
828 
void ComputeMatrix(PROBLEM &pde, MATRIX &matrix)
Definition: integratormixeddims.h:386
const std::map< std::string, const VECTOR * > & GetDomainData() const
Definition: integratormixeddims.h:750
IntegratorMixedDimensions(INTEGRATORDATACONT &idc)
Definition: integratormixeddims.h:135
const INTEGRATORDATACONT & GetIntegratorDataContainer() const
Definition: integratormixeddims.h:803
~IntegratorMixedDimensions()
Definition: integratormixeddims.h:144
const std::map< std::string, const dealii::Vector< SCALAR > * > & GetParamData() const
Definition: integratormixeddims.h:794
void ApplyNewtonBoundaryValues(PROBLEM &pde, VECTOR &u)
Definition: integratormixeddims.h:695
void ComputeNonlinearRhs(PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
Definition: integratormixeddims.h:278
virtual void value(const dealii::Point< dealdim > &p, const unsigned int component, const unsigned int dof_number, dealii::Vector< double > &local_vector) const =0
void ReInit()
Definition: integratormixeddims.h:153
void DeleteDomainData(std::string name)
Definition: integratormixeddims.h:733
void ComputeLocalControlConstraints(PROBLEM &pde, VECTOR &constraints)
Definition: integratormixeddims.h:397
Definition: facedatacontainer.h:48
SCALAR ComputeAlgebraicScalar(PROBLEM &pde)
Definition: integratormixeddims.h:608
void AddPresetRightHandSide(double s, dealii::Vector< SCALAR > &residual) const
Definition: integratormixeddims.h:813
SCALAR ComputeDomainScalar(PROBLEM &pde)
Definition: integratormixeddims.h:409
void DeleteParamData(std::string name)
Definition: integratormixeddims.h:776
void ApplyInitialBoundaryValues(PROBLEM &pde, VECTOR &u)
Definition: integratormixeddims.h:618
Definition: integratormixeddims.h:70
void ApplyTransposedInitialBoundaryValues(PROBLEM &pde, VECTOR &u, SCALAR scale)
Definition: integratormixeddims.h:628
void ComputeNonlinearResidual(PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
Definition: integratormixeddims.h:163
SCALAR ComputeFaceScalar(PROBLEM &pde)
Definition: integratormixeddims.h:545
Definition: transposeddirichletdatainterface.h:37
void AddDomainData(std::string name, const VECTOR *new_data)
Definition: integratormixeddims.h:717
Definition: dopeexception.h:35
SCALAR ComputeBoundaryScalar(PROBLEM &pde)
Definition: integratormixeddims.h:471
SCALAR ComputePointScalar(PROBLEM &pde)
Definition: integratormixeddims.h:451
void AddParamData(std::string name, const dealii::Vector< SCALAR > *new_data)
Definition: integratormixeddims.h:759