DOpE
multimesh_facedatacontainer.h
Go to the documentation of this file.
1 
24 #ifndef MULTIMESH_FACEDATACONTAINER_H_
25 #define MULTIMESH_FACEDATACONTAINER_H_
26 
27 #include "spacetimehandler.h"
28 #include "statespacetimehandler.h"
29 #include "fevalues_wrapper.h"
30 #include "dopeexception.h"
31 
32 #include <deal.II/dofs/dof_handler.h>
33 #include <deal.II/hp/dof_handler.h>
34 
35 using namespace dealii;
36 
37 namespace DOpE
38 {
52  template<template<int, int> class DH, typename VECTOR, int dim>
54  {
55  public:
57  {
58  throw(DOpEException(
59  "Dummy class, this constructor should never get called.",
60  "FaceDataContainer<dealii::DoFHandler , VECTOR, dim>::FaceDataContainer"));
61  }
62  ;
63  };
64 
74  template<typename VECTOR, int dim>
75  class Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>
76  {
77 
78  public:
98  template<template<int, int> class FE, typename SPARSITYPATTERN, int dopedim,
99  int dealdim>
101  const Quadrature<dim - 1>& quad,
102  UpdateFlags update_flags,
103  SpaceTimeHandler<FE, dealii::DoFHandler, SPARSITYPATTERN,
104  VECTOR, dopedim, dealdim>& sth,
105  const typename std::vector<typename dealii::DoFHandler<dim>::cell_iterator>& element,
106  const typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element,
107  const std::map<std::string, const Vector<double>*> &param_values,
108  const std::map<std::string, const VECTOR*> &domain_values, bool/* just for compatibility*/) :
109  param_values_(param_values),
110  domain_values_(domain_values),
111  element_(element),
112  tria_element_(tria_element),
113  state_fe_values_(sth.GetFESystem("state"), quad,
114  update_flags),
115  control_fe_values_(sth.GetFESystem("control"), quad,
116  update_flags)
117  {
118  state_index_ = sth.GetStateIndex();
119  if (state_index_ == 1)
120  control_index_ = 0;
121  else
122  control_index_ = 1;
123  n_q_points_per_element_ = quad.size();
124  n_dofs_per_element_ = element[0]->get_fe().dofs_per_cell;
125  control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
126  state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
127  }
128 
130  {
131  }
132  /*********************************************/
133  /*
134  * This function reinits the FEValues on the actual face. Should
135  * be called prior to any of the get-functions.
136  *
137  * @param face_no The 'local number' (i.e. from the perspective of the actual element) of the
138  * actual face.
139  */
140  inline void
141  ReInit(unsigned int coarse_index,unsigned int fine_index, const FullMatrix<double>& prolongation_matrix, unsigned int face_no);
142 
143  /*********************************************/
147  inline void
149 
150  /*********************************************/
155  inline unsigned int
156  GetNDoFsPerElement() const;
157  inline unsigned int
158  GetNbrNDoFsPerElement() const;
159  inline unsigned int
160  GetNQPoints() const;
161  inline unsigned int
162  GetMaterialId() const;
163  inline unsigned int
164  GetNbrMaterialId() const;
165  inline unsigned int
166  GetNbrMaterialId(unsigned int face) const;
167  inline bool
168  GetIsAtBoundary() const;
169  inline double
170  GetElementDiameter() const;
171  inline unsigned int
172  GetBoundaryIndicator() const;
173  inline const DOpEWrapper::FEFaceValues<dim>&
174  GetFEFaceValuesState() const;
175  inline const DOpEWrapper::FEFaceValues<dim>&
176  GetFEFaceValuesControl() const;
177 
178  /**********************************************/
179  /*
180  * Looks up the given name in parameter_data_ and returns the corresponding value
181  * through 'value'.
182  */
183  void
184  GetParamValues(std::string name, Vector<double>& value) const;
185 
186  /*********************************************/
191  /*
192  * Writes the values of the state variable at the quadrature points into values.
193  */
194  inline void
195  GetFaceValuesState(std::string name, std::vector<double>& values) const;
196  /*********************************************/
197  /*
198  * Same as above for the Vector valued case.
199  */
200  inline void
201  GetFaceValuesState(std::string name, std::vector<Vector<double> >& values) const;
202 
203  /*********************************************/
204  /*
205  * Writes the values of the control variable at the quadrature points into values
206  */
207  inline void
208  GetFaceValuesControl(std::string name, std::vector<double>& values) const;
209  /*********************************************/
210  /*
211  * Same as above for the Vector valued case.
212  */
213  inline void
214  GetFaceValuesControl(std::string name, std::vector<Vector<double> >& values) const;
215  /*********************************************/
216  /*
217  * Writes the values of the state gradient at the quadrature points into values.
218  */
219 
220  template<int targetdim>
221  inline void
222  GetFaceGradsState(std::string name,
223  std::vector<Tensor<1, targetdim> >& values) const;
224 
225  /*********************************************/
226  /*
227  * Same as above for the Vector valued case.
228  */
229  template<int targetdim>
230  inline void
231  GetFaceGradsState(std::string name,
232  std::vector<std::vector<Tensor<1, targetdim> > >& values) const;
233 
234  /*********************************************/
235  /*
236  * Writes the values of the control gradient at the quadrature points into values.
237  */
238  template<int targetdim>
239  inline void
240  GetFaceGradsControl(std::string name,
241  std::vector<Tensor<1, targetdim> >& values) const;
242 
243  /*********************************************/
244  /*
245  * Same as above for the Vector valued case.
246  */
247  template<int targetdim>
248  inline void
249  GetFaceGradsControl(std::string name,
250  std::vector<std::vector<Tensor<1, targetdim> > >& values) const;
251 
252  private:
253  /*
254  * Helper Functions
255  */
256  unsigned int
257  GetStateIndex() const;
258  unsigned int
259  GetControlIndex() const;
260  unsigned int
261  GetCoarseIndex() const { return coarse_index_; }
262  unsigned int
263  GetFineIndex() const { return fine_index_; }
264  /***********************************************************/
268  inline void
269  GetValues(typename dealii::DoFHandler<dim>::cell_iterator element,
270  const FullMatrix<double>& prolongation,
271  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
272  std::vector<double>& values) const;
273  /***********************************************************/
277  inline void
278  GetValues(typename dealii::DoFHandler<dim>::cell_iterator element,
279  const FullMatrix<double>& prolongation,
280  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
281  std::vector<Vector<double> >& values) const;
282  /***********************************************************/
286  template<int targetdim>
287  inline void
288  GetGrads(typename dealii::DoFHandler<dim>::cell_iterator element,
289  const FullMatrix<double>& prolongation,
290  const DOpEWrapper::FEFaceValues<dim>& fe_values,
291  std::string name, std::vector<Tensor<1, targetdim> >& values) const;
292  /***********************************************************/
296  template<int targetdim>
297  inline void
298  GetGrads(typename dealii::DoFHandler<dim>::cell_iterator element,
299  const FullMatrix<double>& prolongation,
300  const DOpEWrapper::FEFaceValues<dim>& fe_values,
301  std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values) const;
302  /***********************************************************/
303  inline const std::map<std::string, const VECTOR*> &
304  GetDomainValues() const
305  {
306  return domain_values_;
307  }
308  /***********************************************************/
309  //"global" member data, part of every instantiation
310  const std::map<std::string, const Vector<double>*> &param_values_;
311  const std::map<std::string, const VECTOR*> &domain_values_;
312  unsigned int state_index_;
313  unsigned int control_index_;
314 
315  const typename std::vector<typename dealii::DoFHandler<dim>::cell_iterator>& element_;
316  const typename std::vector<typename dealii::Triangulation<dim>::cell_iterator>& tria_element_;
317  DOpEWrapper::FEFaceValues<dim> state_fe_values_;
318  DOpEWrapper::FEFaceValues<dim> control_fe_values_;
319 
320  FullMatrix<double> control_prolongation_;
321  FullMatrix<double> state_prolongation_;
322  unsigned int coarse_index_, fine_index_;
323 
324  unsigned int n_q_points_per_element_;
325  unsigned int n_dofs_per_element_;
326 
327  unsigned int face_;
328  };
329 
330 
331  /***********************************************************************/
332  /************************IMPLEMENTATION*for*DoFHandler*********************************/
333  /***********************************************************************/
334 
335  template<typename VECTOR, int dim>
336  void
338  unsigned int coarse_index,
339  unsigned int fine_index,
340  const FullMatrix<double>& prolongation_matrix,
341  unsigned int face_no)
342  {
343  face_ = face_no;
344  coarse_index_ = coarse_index;
345  fine_index_ = fine_index;
346  assert(this->GetControlIndex() < element_.size());
347 
348  if(coarse_index == this->GetStateIndex())
349  {
350  state_prolongation_ = prolongation_matrix;
351  control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
352  }
353  else
354  {
355  if(coarse_index == this->GetControlIndex())
356  {
357  control_prolongation_ = prolongation_matrix;
358  state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
359  }
360  else
361  {
362  control_prolongation_ = IdentityMatrix(element_[this->GetControlIndex()]->get_fe().dofs_per_cell);
363  state_prolongation_ = IdentityMatrix(element_[this->GetStateIndex()]->get_fe().dofs_per_cell);
364  fine_index_ = 0;
365  }
366  }
367 
368  state_fe_values_.reinit(tria_element_[GetFineIndex()], face_no);
369  control_fe_values_.reinit(tria_element_[GetFineIndex()], face_no);
370  }
371  /***********************************************************************/
372  template<typename VECTOR, int dim>
373  unsigned int
375  {
376  return n_dofs_per_element_;
377  }
378 
379  /***********************************************************************/
380 
381  template<typename VECTOR, int dim>
382  unsigned int
384  {
385  throw DOpEException("This function has not been written since we do not know what the right neigbour is!",
386  "Multimesh_FaceDataContainer::GetNbrNDoFsPerElement");
387  }
388 
389  /**********************************************/
390  template<typename VECTOR, int dim>
391  unsigned int
393  {
394  return n_q_points_per_element_;
395  }
396 
397  /**********************************************/
398  template<typename VECTOR, int dim>
399  unsigned int
401  {
402  return tria_element_[GetFineIndex()]->material_id();
403  }
404 
405  /**********************************************/
406  template<typename VECTOR, int dim>
407  unsigned int
409  {
410  return this->GetNbrMaterialId(face_);
411  }
412 
413  template<typename VECTOR, int dim>
414  unsigned int
416  unsigned int /*face*/) const
417  {
418  throw DOpEException("This function has not been written since we do not know what the right neigbour is!",
419  "Multimesh_FaceDataContainer::GetNbrMaterialId");
420  }
421 
422  /**********************************************/
423  template<typename VECTOR, int dim>
424  bool
426  {
427  return tria_element_[GetFineIndex()]->face(face_)->at_boundary();
428  }
429 
430  /**********************************************/
431  template<typename VECTOR, int dim>
432  double
434  {
435  throw DOpEException("This function has not been written since we do not know what the right Diameter!",
436  "Multimesh_FaceDataContainer::GetElementDiameter");
437  }
438 
439  /**********************************************/
440 
441  template<typename VECTOR, int dim>
442  unsigned int
444  {
445  return tria_element_[GetFineIndex()]->face(face_)->boundary_indicator();
446  }
447 
448  /**********************************************/
449  template<typename VECTOR, int dim>
452  {
453  return state_fe_values_;
454  }
455 
456  /**********************************************/
457  template<typename VECTOR, int dim>
460  {
461  return control_fe_values_;
462  }
463 
464  /**********************************************/
465 
466  template<typename VECTOR, int dim>
467  void
469  std::string name, Vector<double>& value) const
470  {
471  typename std::map<std::string, const Vector<double>*>::const_iterator it =
472  param_values_.find(name);
473  if (it == param_values_.end())
474  {
475  throw DOpEException("Did not find " + name,
476  "Multimesh_FaceDataContainer::GetParamValues");
477  }
478  value = *(it->second);
479  }
480 
481  /*********************************************/
482  template<typename VECTOR, int dim>
483  void
485  std::string name, std::vector<double>& values) const
486  {
487  this->GetValues(element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
488  }
489  /*********************************************/
490  template<typename VECTOR, int dim>
491  void
493  std::string name, std::vector<Vector<double> >& values) const
494  {
495  this->GetValues(element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
496 
497  }
498 
499  /*********************************************/
500  template<typename VECTOR, int dim>
501  void
503  std::string name, std::vector<double>& values) const
504  {
505  this->GetValues(element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
506  }
507 
508  /*********************************************/
509  template<typename VECTOR, int dim>
510  void
512  std::string name, std::vector<Vector<double> >& values) const
513  {
514  this->GetValues(element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
515  }
516 
517  /*********************************************/
518  template<typename VECTOR, int dim>
519  template<int targetdim>
520  void
522  std::string name, std::vector<Tensor<1, targetdim> >& values) const
523  {
524  this->GetGrads<targetdim> (element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
525  }
526 
527  /*********************************************/
528  template<typename VECTOR, int dim>
529  template<int targetdim>
530  void
532  std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values) const
533  {
534  this->GetGrads<targetdim> (element_[this->GetStateIndex()],state_prolongation_,this->GetFEFaceValuesState(), name, values);
535  }
536 
537  /***********************************************************************/
538 
539  template<typename VECTOR, int dim>
540  template<int targetdim>
541  void
543  std::string name, std::vector<Tensor<1, targetdim> >& values) const
544  {
545  this->GetGrads<targetdim> (element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
546  }
547  /***********************************************************************/
548 
549  template<typename VECTOR, int dim>
550  template<int targetdim>
551  void
553  std::string name, std::vector<std::vector<Tensor<1, targetdim> > >& values) const
554  {
555  this->GetGrads<targetdim> (element_[this->GetControlIndex()],control_prolongation_,this->GetFEFaceValuesControl(), name, values);
556  }
557 
558  /***********************************************************************/
559 
560  template<typename VECTOR, int dim>
561  unsigned int
563  {
564  return state_index_;
565  }
566 
567  /***********************************************************************/
568 
569  template<typename VECTOR, int dim>
570  unsigned int
571  Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetControlIndex() const
572  {
573  return control_index_;
574  }
575 
576  /***********************************************************************/
577  template<typename VECTOR, int dim>
578  void
579  Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetValues(
580  typename dealii::DoFHandler<dim>::cell_iterator element,
581  const FullMatrix<double>& prolongation,
582  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
583  std::vector<double>& values) const
584  {
585  typename std::map<std::string, const VECTOR*>::const_iterator it =
586  this->GetDomainValues().find(name);
587  if (it == this->GetDomainValues().end())
588  {
589  throw DOpEException("Did not find " + name,
590  "Multimesh_FaceDataContainer::GetValues");
591  }
592  unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
593  //Now we get the values on the real element
594  dealii::Vector<double> dof_values(dofs_per_element);
595  dealii::Vector<double> dof_values_transformed(dofs_per_element);
596  element->get_dof_values (*(it->second), dof_values);
597  //Now compute the real values at the nodal points
598  prolongation.vmult(dof_values_transformed,dof_values);
599 
600  //Copied from deal FEValuesBase<dim,spacedim>::get_function_values
601  // see deal.II/source/fe/fe_values.cc
602  unsigned int n_quadrature_points = GetNQPoints();
603  std::fill_n (values.begin(), n_quadrature_points, 0);
604  for (unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
605  {
606  const double value = dof_values_transformed(shape_func);
607  if (value == 0.)
608  continue;
609 
610  const double *shape_value_ptr = &(fe_values.shape_value(shape_func, 0));
611  for (unsigned int point=0; point<n_quadrature_points; ++point)
612  values[point] += value * *shape_value_ptr++;
613  }
614  }
615 
616  /***********************************************************************/
617  template<typename VECTOR, int dim>
618  void
619  Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetValues(
620  typename dealii::DoFHandler<dim>::cell_iterator element,
621  const FullMatrix<double>& prolongation,
622  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
623  std::vector<Vector<double> >& values) const
624  {
625  typename std::map<std::string, const VECTOR*>::const_iterator it =
626  this->GetDomainValues().find(name);
627  if (it == this->GetDomainValues().end())
628  {
629  throw DOpEException("Did not find " + name,
630  "Multimesh_FaceDataContainer::GetValues");
631  }
632 
633  unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
634  //Now we get the values on the real element
635  dealii::Vector<double> dof_values(dofs_per_element);
636  dealii::Vector<double> dof_values_transformed(dofs_per_element);
637  element->get_dof_values (*(it->second), dof_values);
638  //Now compute the real values at the nodal points
639  prolongation.vmult(dof_values_transformed,dof_values);
640 
641  //Copied from deal FEValuesBase<dim,spacedim>::get_function_values
642  // see deal.II/source/fe/fe_values.cc
643  const unsigned int n_components = element->get_fe().n_components();
644  unsigned int n_quadrature_points = GetNQPoints();
645  for (unsigned i=0;i<values.size();++i)
646  std::fill_n (values[i].begin(), values[i].size(), 0);
647 
648  for (unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
649  {
650  const double value = dof_values_transformed(shape_func);
651  if (value == 0.)
652  continue;
653 
654  if (element->get_fe().is_primitive(shape_func))
655  {
656  const unsigned int comp = element->get_fe().system_to_component_index(shape_func).first;
657  for (unsigned int point=0; point<n_quadrature_points; ++point)
658  values[point](comp) += value * fe_values.shape_value(shape_func,point);
659  }
660  else
661  for (unsigned int c=0; c<n_components; ++c)
662  {
663  if (element->get_fe().get_nonzero_components(shape_func)[c] == false)
664  continue;
665  for (unsigned int point=0; point<n_quadrature_points; ++point)
666  values[point](c) += value * fe_values.shape_value_component(shape_func,point,c);
667  }
668  }
669  }
670 
671  /***********************************************************************/
672 
673  template<typename VECTOR, int dim>
674  template<int targetdim>
675  void
676  Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetGrads(
677  typename dealii::DoFHandler<dim>::cell_iterator element,
678  const FullMatrix<double>& prolongation,
679  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
680  std::vector<Tensor<1, targetdim> >& values) const
681  {
682  typename std::map<std::string, const VECTOR*>::const_iterator it =
683  this->GetDomainValues().find(name);
684  if (it == this->GetDomainValues().end())
685  {
686  throw DOpEException("Did not find " + name,
687  "Multimesh_FaceDataContainerBase::GetGrads");
688  }
689  unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
690  //Now we get the values on the real element
691  dealii::Vector<double> dof_values(dofs_per_element);
692  dealii::Vector<double> dof_values_transformed(dofs_per_element);
693  element->get_dof_values (*(it->second), dof_values);
694  //Now compute the real values at the nodal points
695  prolongation.vmult(dof_values_transformed,dof_values);
696 
697  //Copied from deal FEValuesBase<dim,spacedim>::get_function_gradients
698  unsigned int n_quadrature_points = GetNQPoints();
699  std::fill_n (values.begin(), n_quadrature_points, Tensor<1,targetdim>());
700 
701  for (unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
702  {
703  const double value = dof_values_transformed(shape_func);
704  if (value == 0.)
705  continue;
706 
707  const Tensor<1,targetdim> *shape_gradient_ptr
708  = &(fe_values.shape_grad(shape_func,0));
709  for (unsigned int point=0; point<n_quadrature_points; ++point)
710  values[point] += value * *shape_gradient_ptr++;
711  }
712  }
713 
714  /***********************************************************************/
715 
716  template<typename VECTOR, int dim>
717  template<int targetdim>
718  void
719  Multimesh_FaceDataContainer<dealii::DoFHandler, VECTOR, dim>::GetGrads(
720  typename dealii::DoFHandler<dim>::cell_iterator element,
721  const FullMatrix<double>& prolongation,
722  const DOpEWrapper::FEFaceValues<dim>& fe_values, std::string name,
723  std::vector<std::vector<Tensor<1, targetdim> > >& values) const
724  {
725  typename std::map<std::string, const VECTOR*>::const_iterator it =
726  this->GetDomainValues().find(name);
727  if (it == this->GetDomainValues().end())
728  {
729  throw DOpEException("Did not find " + name,
730  "Multimesh_FaceDataContainerBase::GetGrads");
731  }
732 
733  unsigned int dofs_per_element = element->get_fe().dofs_per_cell;
734  //Now we get the values on the real element
735  dealii::Vector<double> dof_values(dofs_per_element);
736  dealii::Vector<double> dof_values_transformed(dofs_per_element);
737  element->get_dof_values (*(it->second), dof_values);
738  //Now compute the real values at the nodal points
739  prolongation.vmult(dof_values_transformed,dof_values);
740 
741  //Copied from deal FEValuesBase<dim,spacedim>::get_function_gradients
742  const unsigned int n_components = element->get_fe().n_components();
743  unsigned int n_quadrature_points = GetNQPoints();
744  for (unsigned i=0;i<values.size();++i)
745  std::fill_n (values[i].begin(), values[i].size(), Tensor<1,dim>());
746 
747  for (unsigned int shape_func=0; shape_func<dofs_per_element; ++shape_func)
748  {
749  const double value = dof_values_transformed(shape_func);
750  if (value == 0.)
751  continue;
752 
753  if (element->get_fe().is_primitive(shape_func))
754  {
755  const unsigned int comp = element->get_fe().system_to_component_index(shape_func).first;
756  for (unsigned int point=0; point<n_quadrature_points; ++point)
757  values[point][comp] += value * fe_values.shape_grad(shape_func,point);
758  }
759  else
760  for (unsigned int c=0; c<n_components; ++c)
761  {
762  if (element->get_fe().get_nonzero_components(shape_func)[c] == false)
763  continue;
764  for (unsigned int point=0; point<n_quadrature_points; ++point)
765  values[point][c] += value * fe_values.shape_grad_component(shape_func,point,c);
766  }
767  }
768  }
769 
770  /***********************************************************************/
771  /************************END*OF*IMPLEMENTATION**************************/
772  /***********************************************************************/
773 
774 }//end of namespace
775 
776 #endif /* MULTIMESH_FACEDATACONTAINER_H_ */
~Multimesh_FaceDataContainer()
Definition: multimesh_facedatacontainer.h:129
Multimesh_FaceDataContainer(const Quadrature< dim-1 > &quad, UpdateFlags update_flags, SpaceTimeHandler< FE, dealii::DoFHandler, SPARSITYPATTERN, VECTOR, dopedim, dealdim > &sth, const typename std::vector< typename dealii::DoFHandler< dim >::cell_iterator > &element, const typename std::vector< typename dealii::Triangulation< dim >::cell_iterator > &tria_element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, bool)
Definition: multimesh_facedatacontainer.h:100
void ReInitNbr()
Definition: multimesh_facedatacontainer.h:148
Multimesh_FaceDataContainer()
Definition: multimesh_facedatacontainer.h:56
Definition: spacetimehandler.h:71
Definition: multimesh_facedatacontainer.h:53
Definition: fevalues_wrapper.h:88
Definition: dopeexception.h:35