DOpE
functionalinterface.h
Go to the documentation of this file.
1 
24 #ifndef FUNCTIONAL_INTERFACE_H_
25 #define FUNCTIONAL_INTERFACE_H_
26 
27 #include <map>
28 #include <string>
29 
30 #include <deal.II/fe/fe_system.h>
31 #include <deal.II/fe/fe_values.h>
32 #include <deal.II/fe/mapping.h>
33 
34 #include "fevalues_wrapper.h"
35 #include "dofhandler_wrapper.h"
36 #include "elementdatacontainer.h"
37 #include "facedatacontainer.h"
40 
41 namespace DOpE
42 {
48  template<
49  template<template<int, int> class DH, typename VECTOR, int dealdim> class EDC,
50  template<template<int, int> class DH, typename VECTOR, int dealdim> class FDC,
51  template<int, int> class DH, typename VECTOR, int dopedim, int dealdim =
52  dopedim>
54  {
55  public:
57  virtual
59 
66  virtual double
67  ElementValue(const EDC<DH, VECTOR, dealdim>& edc);
68 
79  virtual void
80  ElementValue_U(const EDC<DH, VECTOR, dealdim>& edc,
81  dealii::Vector<double> &local_vector, double scale);
82 
93  virtual void
94  ElementValue_Q(const EDC<DH, VECTOR, dealdim>& edc,
95  dealii::Vector<double> &local_vector, double scale);
96 
107  virtual void
108  ElementValue_UU(const EDC<DH, VECTOR, dealdim>& edc,
109  dealii::Vector<double> &local_vector, double scale);
110 
121  virtual void
122  ElementValue_QU(const EDC<DH, VECTOR, dealdim>& edc,
123  dealii::Vector<double> &local_vector, double scale);
124 
135  virtual void
136  ElementValue_UQ(const EDC<DH, VECTOR, dealdim>& edc,
137  dealii::Vector<double> &local_vector, double scale);
138 
149  virtual void
150  ElementValue_QQ(const EDC<DH, VECTOR, dealdim>& edc,
151  dealii::Vector<double> &local_vector, double scale);
152 
165  virtual double
166  PointValue(
167  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
168  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
169  const std::map<std::string, const dealii::Vector<double>*> &param_values,
170  const std::map<std::string, const VECTOR*> &domain_values);
171 
185  virtual void
186  PointValue_U(
187  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
188  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
189  const std::map<std::string, const dealii::Vector<double>*> &param_values,
190  const std::map<std::string, const VECTOR*> &domain_values,
191  VECTOR& rhs, double scale);
192 
206  virtual void
207  PointValue_Q(
208  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
209  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
210  const std::map<std::string, const dealii::Vector<double>*> &param_values,
211  const std::map<std::string, const VECTOR*> &domain_values,
212  VECTOR& rhs, double scale);
213 
227  virtual void
229  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
230  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
231  const std::map<std::string, const dealii::Vector<double>*> &param_values,
232  const std::map<std::string, const VECTOR*> &domain_values,
233  VECTOR& rhs, double scale);
234 
249  virtual void
251  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
252  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
253  const std::map<std::string, const dealii::Vector<double>*> &param_values,
254  const std::map<std::string, const VECTOR*> &domain_values,
255  VECTOR& rhs, double scale);
256 
271  virtual void
273  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
274  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
275  const std::map<std::string, const dealii::Vector<double>*> &param_values,
276  const std::map<std::string, const VECTOR*> &domain_values,
277  VECTOR& rhs, double scale);
278 
293  virtual void
295  const DOpEWrapper::DoFHandler<dopedim, DH> &control_dof_handler,
296  const DOpEWrapper::DoFHandler<dealdim, DH> &state_dof_handler,
297  const std::map<std::string, const dealii::Vector<double>*> &param_values,
298  const std::map<std::string, const VECTOR*> &domain_values,
299  VECTOR& rhs, double scale);
300 
309  virtual double
310  BoundaryValue(const FDC<DH, VECTOR, dealdim>& fdc);
311 
320  virtual void
321  BoundaryValue_U(const FDC<DH, VECTOR, dealdim>& fdc,
322  dealii::Vector<double> &local_vector, double scale);
323 
332  virtual void
333  BoundaryValue_Q(const FDC<DH, VECTOR, dealdim>& fdc,
334  dealii::Vector<double> &local_vector, double scale);
335 
344  virtual void
345  BoundaryValue_UU(const FDC<DH, VECTOR, dealdim>& fdc,
346  dealii::Vector<double> &local_vector, double scale);
347 
356  virtual void
357  BoundaryValue_QU(const FDC<DH, VECTOR, dealdim>& fdc,
358  dealii::Vector<double> &local_vector, double scale);
359 
368  virtual void
369  BoundaryValue_UQ(const FDC<DH, VECTOR, dealdim>& fdc,
370  dealii::Vector<double> &local_vector, double scale);
371 
380  virtual void
381  BoundaryValue_QQ(const FDC<DH, VECTOR, dealdim>& fdc,
382  dealii::Vector<double> &local_vector, double scale);
383 
391  virtual double
392  FaceValue(const FDC<DH, VECTOR, dealdim>& fdc);
393 
403  virtual void
404  FaceValue_U(const FDC<DH, VECTOR, dealdim>& fdc,
405  dealii::Vector<double> &local_vector, double scale);
406 
416  virtual void
417  FaceValue_Q(const FDC<DH, VECTOR, dealdim>& fdc,
418  dealii::Vector<double> &local_vector, double scale);
419 
429  virtual void
430  FaceValue_UU(const FDC<DH, VECTOR, dealdim>& fdc,
431  dealii::Vector<double> &local_vector, double scale);
432 
442  virtual void
443  FaceValue_QU(const FDC<DH, VECTOR, dealdim>& fdc,
444  dealii::Vector<double> &local_vector, double scale);
445 
455  virtual void
456  FaceValue_UQ(const FDC<DH, VECTOR, dealdim>& fdc,
457  dealii::Vector<double> &local_vector, double scale);
458 
468  virtual void
469  FaceValue_QQ(const FDC<DH, VECTOR, dealdim>& fdc,
470  dealii::Vector<double> &local_vector, double scale);
471 
475  virtual double
477  const std::map<std::string, const dealii::Vector<double>*> &param_values,
478  const std::map<std::string, const VECTOR*> &domain_values);
482  virtual void
483  AlgebraicGradient_Q(VECTOR& gradient,
484  const std::map<std::string, const dealii::Vector<double>*> &param_values,
485  const std::map<std::string, const VECTOR*> &domain_values);
486 
498  virtual std::string
499  GetType() const;
505  virtual std::string
506  GetName() const;
507 
516  virtual bool
517  NeedTime() const
518  {
519  return true;
520  }
521 
528  virtual void
529  SetTime(double /*t*/) const
530  {
531  }
532 
537  virtual dealii::UpdateFlags
538  GetUpdateFlags() const;
539 
544  virtual dealii::UpdateFlags
545  GetFaceUpdateFlags() const;
546 
554  virtual bool
555  HasFaces() const;
556 
564  virtual bool
565  HasPoints() const;
566 
567  void
568  SetProblemType(std::string type);
569  protected:
570  std::string
571  GetProblemType() const;
572  private:
573  std::string problem_type_;
574 
575  };
576 }
577 
578 #endif
virtual double AlgebraicValue(const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: functionalinterface.cc:515
virtual void ElementValue_U(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:198
virtual void BoundaryValue_Q(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:440
virtual void BoundaryValue_QQ(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:500
virtual void FaceValue_QU(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:367
virtual bool NeedTime() const
Definition: functionalinterface.h:517
virtual void PointValue_QU(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:146
virtual void ElementValue_UQ(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:258
virtual void ElementValue_Q(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:213
virtual void ElementValue_QQ(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:273
std::string GetProblemType() const
Definition: functionalinterface.cc:618
virtual void FaceValue_U(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:322
virtual void BoundaryValue_UQ(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:485
virtual void SetTime(double) const
Definition: functionalinterface.h:529
virtual void FaceValue_QQ(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:397
virtual void BoundaryValue_UU(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:455
virtual void BoundaryValue_QU(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:470
virtual dealii::UpdateFlags GetUpdateFlags() const
Definition: functionalinterface.cc:530
virtual std::string GetType() const
Definition: functionalinterface.cc:288
virtual double BoundaryValue(const FDC< DH, VECTOR, dealdim > &fdc)
Definition: functionalinterface.cc:412
virtual void PointValue_Q(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:110
virtual void FaceValue_UU(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:352
virtual void FaceValue_UQ(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:382
virtual void PointValue_U(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:92
virtual dealii::UpdateFlags GetFaceUpdateFlags() const
Definition: functionalinterface.cc:542
Definition: functionalinterface.h:53
virtual void PointValue_UU(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:128
virtual std::string GetName() const
Definition: functionalinterface.cc:299
virtual void AlgebraicGradient_Q(VECTOR &gradient, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: functionalinterface.cc:554
virtual double FaceValue(const FDC< DH, VECTOR, dealdim > &fdc)
Definition: functionalinterface.cc:310
virtual void PointValue_UQ(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:163
void SetProblemType(std::string type)
Definition: functionalinterface.cc:606
virtual double ElementValue(const EDC< DH, VECTOR, dealdim > &edc)
Definition: functionalinterface.cc:63
virtual void FaceValue_Q(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:337
virtual bool HasPoints() const
Definition: functionalinterface.cc:589
virtual ~FunctionalInterface()
Definition: functionalinterface.cc:51
virtual void ElementValue_QU(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:243
FunctionalInterface()
Definition: functionalinterface.cc:40
virtual void BoundaryValue_U(const FDC< DH, VECTOR, dealdim > &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:425
virtual void PointValue_QQ(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs, double scale)
Definition: functionalinterface.cc:180
virtual void ElementValue_UU(const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
Definition: functionalinterface.cc:228
virtual double PointValue(const DOpEWrapper::DoFHandler< dopedim, DH > &control_dof_handler, const DOpEWrapper::DoFHandler< dealdim, DH > &state_dof_handler, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values)
Definition: functionalinterface.cc:76
virtual bool HasFaces() const
Definition: functionalinterface.cc:570