DOpE
elementdatacontainer.h
Go to the documentation of this file.
1 
24 #ifndef ELEMENTDATACONTAINER_H_
25 #define ELEMENTDATACONTAINER_H_
26 
27 #include "spacetimehandler.h"
28 #include "statespacetimehandler.h"
29 #include "fevalues_wrapper.h"
30 #include "dopeexception.h"
32 
33 #include <dofs/dof_handler.h>
34 #include <multigrid/mg_dof_handler.h>
35 #include <hp/dof_handler.h>
36 
37 using namespace dealii;
38 
39 namespace DOpE
40 {
55  template<template<int, int> class DH, typename VECTOR, int dim>
57  VECTOR, dim>
58  {
59  public:
61  {
62  throw(DOpE::DOpEException(
63  "Dummy class, this constructor should never get called.",
64  "ElementDataContainer<dealii::DoFHandler<dim> , VECTOR, dim>::ElementDataContainer"));
65  }
66  ;
67  };
68 
78  template<typename VECTOR, int dim>
79  class ElementDataContainer<dealii::DoFHandler, VECTOR, dim> : public edcinternal::ElementDataContainerInternal<
80  VECTOR, dim>
81  {
82 
83  public:
104  template<template<int, int> class FE, typename SPARSITYPATTERN, int dopedim, int dealdim>
105  ElementDataContainer(const Quadrature<dim>& quad,
106  UpdateFlags update_flags,
107  SpaceTimeHandler<FE, dealii::DoFHandler, SPARSITYPATTERN, VECTOR,
108  dopedim, dealdim>& sth,
109  const std::vector<
110  typename dealii::DoFHandler<dim>::active_cell_iterator>& element,
111  const std::map<std::string, const Vector<double>*> &param_values,
112  const std::map<std::string, const VECTOR*> &domain_values) :
113  edcinternal::ElementDataContainerInternal<VECTOR, dim>(param_values,
114  domain_values), _element(element), _state_fe_values(
115  sth.GetMapping(), (sth.GetFESystem("state")), quad,
116  update_flags), _control_fe_values(sth.GetMapping(),
117  (sth.GetFESystem("control")), quad, update_flags)
118  {
119  _state_index = sth.GetStateIndex();
120  if (_state_index == 1)
121  _control_index = 0;
122  else
123  _control_index = 1;
124  _n_q_points_per_element = quad.size();
125  _n_dofs_per_element = element[0]->get_fe().dofs_per_cell;
126  }
127 
146  template<template<int, int> class FE, typename SPARSITYPATTERN>
147  ElementDataContainer(const Quadrature<dim>& quad,
148  UpdateFlags update_flags,
149  StateSpaceTimeHandler<FE, dealii::DoFHandler, SPARSITYPATTERN,
150  VECTOR, dim>& sth,
151  const std::vector<
152  typename dealii::DoFHandler<dim>::active_cell_iterator>& element,
153  const std::map<std::string, const Vector<double>*> &param_values,
154  const std::map<std::string, const VECTOR*> &domain_values) :
155  edcinternal::ElementDataContainerInternal<VECTOR, dim>(param_values,
156  domain_values), _element(element), _state_fe_values(
157  sth.GetMapping(), (sth.GetFESystem("state")), quad,
158  update_flags), _control_fe_values(sth.GetMapping(),
159  (sth.GetFESystem("state")), quad, update_flags)
160  {
161  _state_index = sth.GetStateIndex();
162  _control_index = element.size(); //Make sure they are never used ...
163  _n_q_points_per_element = quad.size();
164  _n_dofs_per_element = element[0]->get_fe().dofs_per_cell;
165  }
167  {
168  }
169  /*********************************************/
170  /*
171  * This function reinits the FEValues on the actual element. Should
172  * be called prior to any of the get-functions.
173  */
174  inline void
175  ReInit();
176 
177  /*********************************************/
182  inline unsigned int
183  GetNDoFsPerElement() const;
184  inline unsigned int
185  GetNQPoints() const;
186  inline unsigned int
187  GetMaterialId() const;
188  inline unsigned int
189  GetNbrMaterialId(unsigned int face) const;
190  inline unsigned int
191  GetFaceBoundaryIndicator(unsigned int face) const;
192  inline bool
193  GetIsAtBoundary() const;
194  inline double
195  GetElementDiameter() const;
196  inline Point<dim>
197  GetCenter() const;
198  inline const DOpEWrapper::FEValues<dim>&
199  GetFEValuesState() const;
200  inline const DOpEWrapper::FEValues<dim>&
201  GetFEValuesControl() const;
202  private:
203  /*
204  * Helper Functions
205  */
206  unsigned int
207  GetStateIndex() const;
208  unsigned int
209  GetControlIndex() const;
210 
211  /***********************************************************/
212  //"global" member data, part of every instantiation
213  unsigned int _state_index;
214  unsigned int _control_index;
215 
216  const std::vector<typename dealii::DoFHandler<dim>::active_cell_iterator> & _element;
217  DOpEWrapper::FEValues<dim> _state_fe_values;
218  DOpEWrapper::FEValues<dim> _control_fe_values;
219 
220  unsigned int _n_q_points_per_element;
221  unsigned int _n_dofs_per_element;
222  };
223 
224  /*****************************************************************/
225  /* MGDofHandler */
226 
227 
237  template<typename VECTOR, int dim>
238  class ElementDataContainer<dealii::MGDoFHandler, VECTOR, dim> : public edcinternal::ElementDataContainerInternal<
239  VECTOR, dim>
240  {
241 
242  public:
263  template<template<int, int> class FE, typename SPARSITYPATTERN, int dopedim, int dealdim>
265  const Quadrature<dim>& quad,
266  UpdateFlags update_flags,
267  SpaceTimeHandler<FE, dealii::MGDoFHandler, SPARSITYPATTERN,
268  VECTOR, dopedim, dealdim>& sth
269  ,
270  const std::vector<
271  typename dealii::MGDoFHandler<dim>::active_cell_iterator>& element,
272  const std::map<std::string, const Vector<double>*> &param_values,
273  const std::map<std::string, const VECTOR*> &domain_values)
274  : edcinternal::ElementDataContainerInternal<VECTOR, dim>(
275  param_values, domain_values), _element(element), _state_fe_values(
276  sth.GetMapping(), (sth.GetFESystem("state")), quad,
277  update_flags), _control_fe_values(sth.GetMapping(),
278  (sth.GetFESystem("control")), quad, update_flags)
279  {
280  _state_index = sth.GetStateIndex();
281  if (_state_index == 1)
282  _control_index = 0;
283  else
284  _control_index = 1;
285  _n_q_points_per_element = quad.size();
286  _n_dofs_per_element = element[0]->get_fe().dofs_per_cell;
287  }
288 
289 
290 
291 
292 
311  template<template<int, int> class FE, typename SPARSITYPATTERN>
312  ElementDataContainer(const Quadrature<dim>& quad,
313  UpdateFlags update_flags,
314  StateSpaceTimeHandler<FE, dealii::MGDoFHandler,
315  SPARSITYPATTERN, VECTOR, dim>& sth,
316  const std::vector<
317  typename dealii::MGDoFHandler<dim>::active_cell_iterator>& element,
318  const std::map<std::string, const Vector<double>*> &param_values,
319  const std::map<std::string, const VECTOR*> &domain_values)
320  : edcinternal::ElementDataContainerInternal<VECTOR, dim>(
321  param_values, domain_values), _element(element), _state_fe_values(
322  sth.GetMapping(), (sth.GetFESystem("state")), quad,
323  update_flags), _control_fe_values(sth.GetMapping(),
324  (sth.GetFESystem("state")), quad, update_flags)
325  {
326  _state_index = sth.GetStateIndex();
327  _control_index = element.size(); //Make sure they are never used ...
328  _n_q_points_per_element = quad.size();
329  _n_dofs_per_element = element[0]->get_fe().dofs_per_cell;
330  }
331 
332 
334  {
335  }
336  /*********************************************/
337  /*
338  * This function reinits the FEValues on the actual element. Should
339  * be called prior to any of the get-functions.
340  */
341  inline void
342  ReInit();
343 
344  /*********************************************/
349  inline unsigned int
350  GetNDoFsPerElement() const;
351  inline unsigned int
352  GetNQPoints() const;
353  inline unsigned int
354  GetMaterialId() const;
355  inline unsigned int
356  GetNbrMaterialId(unsigned int face) const;
357  inline unsigned int
358  GetFaceBoundaryIndicator(unsigned int face) const;
359  inline bool
360  GetIsAtBoundary() const;
361  inline double
362  GetElementDiameter() const;
363  inline Point<dim> GetCenter() const;
364  inline const DOpEWrapper::FEValues<dim>&
365  GetFEValuesState() const;
366  inline const DOpEWrapper::FEValues<dim>&
367  GetFEValuesControl() const;
368  private:
369  /*
370  * Helper Functions
371  */
372  unsigned int
373  GetStateIndex() const;
374  unsigned int
375  GetControlIndex() const;
376 
377  /***********************************************************/
378  //"global" member data, part of every instantiation
379  unsigned int _state_index;
380  unsigned int _control_index;
381 
382  const std::vector<typename dealii::MGDoFHandler<dim>::active_cell_iterator> & _element;
383  DOpEWrapper::FEValues<dim> _state_fe_values;
384  DOpEWrapper::FEValues<dim> _control_fe_values;
385 
386  unsigned int _n_q_points_per_element;
387  unsigned int _n_dofs_per_element;
388  };
389 
390  /* end MGDofHandler */
391  /*************************************************************************/
392 
393 
394  template<typename VECTOR, int dim>
395  class ElementDataContainer<dealii::hp::DoFHandler, VECTOR, dim> : public edcinternal::ElementDataContainerInternal<
396  VECTOR, dim>
397  {
398 
399  public:
420  template<template<int, int> class FE, typename SPARSITYPATTERN, int dopedim, int dealdim>
421  ElementDataContainer(const hp::QCollection<dim>& q_collection,
422  UpdateFlags update_flags,
423  SpaceTimeHandler<FE, dealii::hp::DoFHandler, SPARSITYPATTERN,
424  VECTOR, dopedim, dealdim>& sth,
425  const std::vector<
427  const std::map<std::string, const Vector<double>*> &param_values,
428  const std::map<std::string, const VECTOR*> &domain_values) :
429  edcinternal::ElementDataContainerInternal<VECTOR, dim>(param_values,
430  domain_values), _element(element), _state_hp_fe_values(
431  sth.GetMapping(), (sth.GetFESystem("state")), q_collection,
432  update_flags), _control_hp_fe_values(sth.GetMapping(),
433  (sth.GetFESystem("control")), q_collection, update_flags), _q_collection(
434  q_collection)
435  {
436  _state_index = sth.GetStateIndex();
437  if (_state_index == 1)
438  _control_index = 0;
439  else
440  _control_index = 1;
441  }
442 
460  template<template<int, int> class FE, typename SPARSITYPATTERN>
461  ElementDataContainer(const hp::QCollection<dim>& q_collection,
462  UpdateFlags update_flags,
463  StateSpaceTimeHandler<FE, dealii::hp::DoFHandler, SPARSITYPATTERN,
464  VECTOR, dim>& sth,
465  const std::vector<
467  const std::map<std::string, const Vector<double>*> &param_values,
468  const std::map<std::string, const VECTOR*> &domain_values) :
469  edcinternal::ElementDataContainerInternal<VECTOR, dim>(param_values,
470  domain_values), _element(element), _state_hp_fe_values(
471  sth.GetMapping(), (sth.GetFESystem("state")), q_collection,
472  update_flags), _control_hp_fe_values(sth.GetMapping(),
473  (sth.GetFESystem("state")), q_collection, update_flags), _q_collection(
474  q_collection)
475  {
476  _state_index = sth.GetStateIndex();
477  _control_index = element.size(); //Make sure they are never used ...
478  }
480  {
481  }
482 
483  /*****************************************************************/
484  /*
485  * This function reinits the hp::FEValues on the actual element. Should
486  * be called prior to any of the get-functions.
487  */
488  inline void
489  ReInit();
490  /*********************************************/
496  inline unsigned int
497  GetNDoFsPerElement() const;
498  inline unsigned int
499  GetNQPoints() const;
500  inline unsigned int
501  GetMaterialId() const;
502  inline unsigned int
503  GetNbrMaterialId(unsigned int face) const;
504  inline unsigned int
505  GetFaceBoundaryIndicator(unsigned int face) const;
506  inline bool
507  GetIsAtBoundary() const;
508  inline double
509  GetElementDiameter() const;
510  inline Point<dim>
511  GetCenter() const;
512 
513  inline const DOpEWrapper::FEValues<dim>&
514  GetFEValuesState() const;
515  inline const DOpEWrapper::FEValues<dim>&
516  GetFEValuesControl() const;
517 
518  private:
519  unsigned int
520  GetStateIndex() const;
521  unsigned int
522  GetControlIndex() const;
523  const std::map<std::string, const VECTOR*> &
524  GetDomainValues() const;
525  /***********************************************************/
530  inline void
531  GetValues(const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
532  std::vector<double>& values) const;
533  /***********************************************************/
537  inline void
538  GetValues(const DOpEWrapper::FEValues<dim>& fe_values, std::string name,
539  std::vector<Vector<double> >& values) const;
540 
541  /***********************************************************/
545  template<int targetdim>
546  inline void
547  GetGrads(const DOpEWrapper::FEValues<dim>& fe_values,
548  std::string name,
549  std::vector<Tensor<1, targetdim> >& values) const;
550  /***********************************************************/
554  template<int targetdim>
555  inline void
556  GetGrads(const DOpEWrapper::FEValues<dim>& fe_values,
557  std::string name,
558  std::vector<std::vector<Tensor<1, targetdim> > >& values) const;
559 
560  /***********************************************************/
561  //"global" member data, part of every instantiation
562  unsigned int _state_index;
563  unsigned int _control_index;
564 
565  const std::vector<
566  typename dealii::hp::DoFHandler<dim>::active_cell_iterator> & _element;
567  DOpEWrapper::HpFEValues<dim> _state_hp_fe_values;
568  DOpEWrapper::HpFEValues<dim> _control_hp_fe_values;
569 
570  const hp::QCollection<dim>& _q_collection;
571  };
572 
573  /***********************************************************************/
574  /************************IMPLEMENTATION for DoFHandler*********************************/
575  /***********************************************************************/
576 
577  template<typename VECTOR, int dim>
578  void
580  {
581  _state_fe_values.reinit(_element[this->GetStateIndex()]);
582  //Make sure that the Control must be initialized.
583  if (this->GetControlIndex() < _element.size())
584  _control_fe_values.reinit(_element[this->GetControlIndex()]);
585  }
586 
587  /***********************************************************************/
588  template<typename VECTOR, int dim>
589  unsigned int
591  {
592  return _n_dofs_per_element;
593  }
594 
595  /**********************************************/
596  template<typename VECTOR, int dim>
597  unsigned int
599  {
600  return _n_q_points_per_element;
601  }
602 
603  /**********************************************/
604  template<typename VECTOR, int dim>
605  unsigned int
607  {
608  return _element[0]->material_id();
609  }
610 
611  /**********************************************/
612  template<typename VECTOR, int dim>
613  unsigned int
615  unsigned int face) const
616  {
617  if (_element[0]->neighbor_index(face) != -1)
618  return _element[0]->neighbor(face)->material_id();
619  else
620  throw DOpEException("There is no neighbor with number " + face,
621  "ElementDataContainer::GetNbrMaterialId");
622  }
623 
624  /**********************************************/
625  template<typename VECTOR, int dim>
626  unsigned int
628  unsigned int face) const
629  {
630  return _element[0]->face(face)->boundary_indicator();
631  }
632 
633  /**********************************************/
634  template<typename VECTOR, int dim>
635  bool
637  {
638  return _element[0]->at_boundary();
639  }
640  /**********************************************/
641  template<typename VECTOR, int dim>
642  double
644  {
645  return _element[0]->diameter();
646  }
647  /**********************************************/
648  template<typename VECTOR, int dim>
649  Point<dim>
651  {
652  return _element[0]->center();
653  }
654 
655  /**********************************************/
656  template<typename VECTOR, int dim>
659  {
660  return _state_fe_values;
661  }
662 
663  /**********************************************/
664  template<typename VECTOR, int dim>
667  {
668  return _control_fe_values;
669  }
670 
671  /***********************************************************************/
672 
673  template<typename VECTOR, int dim>
674  unsigned int
676  {
677  return _state_index;
678  }
679 
680  /***********************************************************************/
681 
682  template<typename VECTOR, int dim>
683  unsigned int
684  ElementDataContainer<dealii::DoFHandler, VECTOR, dim>::GetControlIndex() const
685  {
686  return _control_index;
687  }
688 
689  /***********************************************************************/
690  /************************END*OF*IMPLEMENTATION**************************/
691  /***********************************************************************/
692  /***********************************************************************/
693  /************************IMPLEMENTATION*********************************/
694  /***********************************************************************/
695 
696 
697 
698  /***********************************************************************/
699  /************************IMPLEMENTATION for MGDoFHandler*********************************/
700  /***********************************************************************/
701 
702  template<typename VECTOR, int dim>
703  void
705  {
706  _state_fe_values.reinit(_element[this->GetStateIndex()]);
707  //Make sure that the Control must be initialized.
708  if (this->GetControlIndex() < _element.size())
709  _control_fe_values.reinit(_element[this->GetControlIndex()]);
710  }
711 
712  /***********************************************************************/
713  template<typename VECTOR, int dim>
714  unsigned int
716  {
717  return _n_dofs_per_element;
718  }
719 
720  /**********************************************/
721  template<typename VECTOR, int dim>
722  unsigned int
724  {
725  return _n_q_points_per_element;
726  }
727 
728  /**********************************************/
729  template<typename VECTOR, int dim>
730  unsigned int
732  {
733  return _element[0]->material_id();
734  }
735 
736  /**********************************************/
737  template<typename VECTOR, int dim>
738  unsigned int
740  unsigned int face) const
741  {
742  if (_element[0]->neighbor_index(face) != -1)
743  return _element[0]->neighbor(face)->material_id();
744  else
745  throw DOpEException("There is no neighbor with number " + face,
746  "ElementDataContainer::GetNbrMaterialId");
747  }
748 
749  /**********************************************/
750  template<typename VECTOR, int dim>
751  unsigned int
753  unsigned int face) const
754  {
755  return _element[0]->face(face)->boundary_indicator();
756  }
757 
758  /**********************************************/
759  template<typename VECTOR, int dim>
760  bool
762  {
763  return _element[0]->at_boundary();
764  }
765  /**********************************************/
766  template<typename VECTOR, int dim>
767  double
769  {
770  return _element[0]->diameter();
771  }
772  /**********************************************/
773  template<typename VECTOR, int dim>
774  Point<dim>
776  {
777  return _element[0]->center();
778  }
779 
780  /**********************************************/
781  template<typename VECTOR, int dim>
784  {
785  return _state_fe_values;
786  }
787 
788  /**********************************************/
789  template<typename VECTOR, int dim>
792  {
793  return _control_fe_values;
794  }
795 
796  /***********************************************************************/
797 
798  template<typename VECTOR, int dim>
799  unsigned int
801  {
802  return _state_index;
803  }
804 
805  /***********************************************************************/
806 
807  template<typename VECTOR, int dim>
808  unsigned int
809  ElementDataContainer<dealii::MGDoFHandler, VECTOR, dim>::GetControlIndex() const
810  {
811  return _control_index;
812  }
813 
814  /***********************************************************************/
815  /************************END*OF*IMPLEMENTATION**************************/
816  /***********************************************************************/
817  /***********************************************************************/
818  /************************IMPLEMENTATION*********************************/
819  /***********************************************************************/
820 
821 
822 
823 
824  template<typename VECTOR, int dim>
825  void
827  {
828  _state_hp_fe_values.reinit(_element[this->GetStateIndex()]);
829  //Make sure that the Control must be initialized.
830  if (this->GetControlIndex() < _element.size())
831  _control_hp_fe_values.reinit(_element[this->GetControlIndex()]);
832  }
833  /*********************************************/
834 
835  template<typename VECTOR, int dim>
836  unsigned int
838  {
839  return _element[0]->get_fe().dofs_per_cell;
840  }
841  /*********************************************/
842  template<typename VECTOR, int dim>
843  unsigned int
845  {
846  return (_q_collection[_element[0]->active_fe_index()]).size();
847  }
848  /*********************************************/
849 
850  template<typename VECTOR, int dim>
851  unsigned int
853  {
854  return _element[0]->material_id();
855  }
856  /*********************************************/
857 
858  template<typename VECTOR, int dim>
859  unsigned int
861  unsigned int face) const
862  {
863  if (_element[0]->neighbor_index(face) != -1)
864  return _element[0]->neighbor(face)->material_id();
865  else
866  throw DOpEException("There is no neighbor with number" + face,
867  "HpElementDataContainer::GetNbrMaterialId");
868  }
869 
870  /**********************************************/
871  template<typename VECTOR, int dim>
872  unsigned int
874  unsigned int face) const
875  {
876  return _element[0]->face(face)->boundary_indicator();
877  }
878 
879  /*********************************************/
880 
881  template<typename VECTOR, int dim>
882  bool
884  {
885  return _element[0]->at_boundary();
886  }
887 
888  /*********************************************/
889 
890  template<typename VECTOR, int dim>
891  double
893  {
894  return _element[0]->diameter();
895  }
896 
897  /**********************************************/
898  template<typename VECTOR, int dim>
899  Point<dim>
901  {
902  return _element[0]->center();
903  }
904 
905  /*********************************************/
906  template<typename VECTOR, int dim>
909  {
910  return static_cast<const DOpEWrapper::FEValues<dim>&>(_state_hp_fe_values.get_present_fe_values());
911  }
912  /*********************************************/
913  template<typename VECTOR, int dim>
916  {
917  return static_cast<const DOpEWrapper::FEValues<dim>&>(_control_hp_fe_values.get_present_fe_values());
918  }
919  /*********************************************/
920 
921  template<typename VECTOR, int dim>
922  unsigned int
924  {
925  return _state_index;
926  }
927  /*********************************************/
928  template<typename VECTOR, int dim>
929  unsigned int
930  ElementDataContainer<dealii::hp::DoFHandler, VECTOR, dim>::GetControlIndex() const
931  {
932  return _control_index;
933  }
934 } //end of namespace
935 
936 #endif /* WORKINGTITLE_H_ */
ElementDataContainer(const Quadrature< dim > &quad, UpdateFlags update_flags, StateSpaceTimeHandler< FE, dealii::MGDoFHandler, SPARSITYPATTERN, VECTOR, dim > &sth, const std::vector< typename dealii::MGDoFHandler< dim >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:312
ElementDataContainer()
Definition: elementdatacontainer.h:60
ElementDataContainer(const hp::QCollection< dim > &q_collection, UpdateFlags update_flags, SpaceTimeHandler< FE, dealii::hp::DoFHandler, SPARSITYPATTERN, VECTOR, dopedim, dealdim > &sth, const std::vector< typename DOpEWrapper::DoFHandler< dim, dealii::hp::DoFHandler >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:421
Definition: elementdatacontainer.h:56
~ElementDataContainer()
Definition: elementdatacontainer.h:166
ElementDataContainer(const Quadrature< dim > &quad, UpdateFlags update_flags, StateSpaceTimeHandler< FE, dealii::DoFHandler, SPARSITYPATTERN, VECTOR, dim > &sth, const std::vector< typename dealii::DoFHandler< dim >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:147
ElementDataContainer(const Quadrature< dim > &quad, UpdateFlags update_flags, SpaceTimeHandler< FE, dealii::DoFHandler, SPARSITYPATTERN, VECTOR, dopedim, dealdim > &sth, const std::vector< typename dealii::DoFHandler< dim >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:105
~ElementDataContainer()
Definition: elementdatacontainer.h:333
Definition: elementdatacontainer_internal.h:46
Definition: dofhandler_wrapper.h:51
Definition: statespacetimehandler.h:59
~ElementDataContainer()
Definition: elementdatacontainer.h:479
Definition: spacetimehandler.h:71
ElementDataContainer(const Quadrature< dim > &quad, UpdateFlags update_flags, SpaceTimeHandler< FE, dealii::MGDoFHandler, SPARSITYPATTERN, VECTOR, dopedim, dealdim > &sth, const std::vector< typename dealii::MGDoFHandler< dim >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:264
ElementDataContainer(const hp::QCollection< dim > &q_collection, UpdateFlags update_flags, StateSpaceTimeHandler< FE, dealii::hp::DoFHandler, SPARSITYPATTERN, VECTOR, dim > &sth, const std::vector< typename DOpEWrapper::DoFHandler< dim, dealii::hp::DoFHandler >::active_cell_iterator > &element, const std::map< std::string, const Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: elementdatacontainer.h:461
Definition: fevalues_wrapper.h:43
Definition: dopeexception.h:35
Definition: fevalues_wrapper.h:186