DOpE
stateproblem.h
Go to the documentation of this file.
1 
24 #ifndef STATE_PROBLEM_H_
25 #define STATE_PROBLEM_H_
26 
27 #include "spacetimehandler.h"
28 
29 using namespace dealii;
30 
31 namespace DOpE
32 {
46  template<typename OPTPROBLEM, typename PDE, typename DD,
47  typename SPARSITYPATTERN, typename VECTOR, int dim>
49  {
50  public:
51  StateProblem(OPTPROBLEM& OP, PDE& pde) :
52  pde_(pde), opt_problem_(OP)
53  {
54  dirichlet_colors_ = opt_problem_.dirichlet_colors_;
55  dirichlet_comps_ = opt_problem_.dirichlet_comps_;
56  primal_dirichlet_values_ = opt_problem_.primal_dirichlet_values_;
57  state_boundary_equation_colors_
58  = opt_problem_.state_boundary_equation_colors_;
59  interval_length_ = 1.;
60  }
61 
62  std::string
63  GetName() const
64  {
65  return "StateProblem";
66  }
67  std::string
68  GetType() const
69  {
70  return "state";
71  }
72 
73  /******************************************************/
74  /****For the initial values ***************/
79  template<typename EDC>
80  void Init_ElementEquation(const EDC& edc,
81  dealii::Vector<double> &local_vector, double scale,
82  double scale_ico)
83  {
84  pde_.Init_ElementEquation(edc, local_vector, scale, scale_ico);
85  }
86 
91  template<typename EDC>
92  void
93  Init_ElementRhs(const EDC& edc,
94  dealii::Vector<double> &local_vector, double scale)
95  {
96  pde_.Init_ElementRhs(& GetInitialValues(), edc, local_vector, scale);
97  }
98 
103  void
105  const std::map<std::string, const dealii::Vector<double>*> &/*param_values*/,
106  const std::map<std::string, const VECTOR*> &/*domain_values*/,
107  VECTOR& /*rhs_vector*/, double /*scale=1.*/)
108  {
109  //Note if this is implemented one needs to update Init_PointRhs in the
110  // OptProblem container in the tangent case.
111  }
112 
117  template<typename EDC>
118  void Init_ElementMatrix(const EDC& edc,
119  dealii::FullMatrix<double> &local_entry_matrix, double scale,
120  double scale_ico)
121  {
122  pde_.Init_ElementMatrix(edc, local_entry_matrix, scale, scale_ico);
123  }
124 
125  /******************************************************/
126  /* Functions as in OptProblem */
131  template<typename EDC>
132  inline void
133  ElementEquation(const EDC& edc,
134  dealii::Vector<double> &local_vector, double scale,
135  double scale_ico);
136 
141  template<typename EDC>
142  inline void
143  ElementTimeEquation(const EDC& edc,
144  dealii::Vector<double> &local_vector, double scale = 1.);
145 
150  template<typename EDC>
151  inline void
152  ElementTimeEquationExplicit(const EDC& edc,
153  dealii::Vector<double> &local_vector, double scale = 1.);
154 
159  template<typename EDC>
160  inline void
161  ElementRhs(const EDC& edc,
162  dealii::Vector<double> &local_vector, double scale = 1.);
163 
168  void
169  PointRhs(
170  const std::map<std::string, const dealii::Vector<double>*> &param_values,
171  const std::map<std::string, const VECTOR*> &domain_values,
172  VECTOR& rhs_vector, double scale);
173 
178  template<typename EDC>
179  inline void
180  ElementMatrix(const EDC& edc,
181  dealii::FullMatrix<double> &local_entry_matrix, double scale = 1.,
182  double scale_ico = 1.);
183 
188  template<typename EDC>
189  inline void
190  ElementTimeMatrix(const EDC& edc,
191  dealii::FullMatrix<double> &local_entry_matrix);
192 
197  template<typename EDC>
198  inline void
199  ElementTimeMatrixExplicit(const EDC& edc,
200  dealii::FullMatrix<double> &local_entry_matrix);
201 
206  template<typename FDC>
207  inline void
208  FaceEquation(const FDC& fdc,
209  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
210 
215  template<typename FDC>
216  inline void
217  InterfaceEquation(const FDC& fdc,
218  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
219 
224  template<typename FDC>
225  inline void
226  FaceRhs(const FDC& fdc,
227  dealii::Vector<double> &local_vector, double scale = 1.);
228 
233  template<typename FDC>
234  inline void
235  FaceMatrix(const FDC& fdc,
236  dealii::FullMatrix<double> &local_entry_matrix, double scale = 1.,double scale_ico = 1.);
237 
238  template<typename FDC>
239  inline void
240  InterfaceMatrix(const FDC& fdc,
241  dealii::FullMatrix<double> &local_entry_matrix, double scale = 1.,double scale_ico = 1.);
242 
247  template<typename FDC>
248  inline void
249  BoundaryEquation(const FDC& fdc,
250  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
251 
256  template<typename FDC>
257  inline void
258  BoundaryRhs(const FDC& fdc,
259  dealii::Vector<double> &local_vector, double scale = 1.);
260 
265  template<typename FDC>
266  inline void
267  BoundaryMatrix(const FDC& fdc,
268  dealii::FullMatrix<double> &local_matrix, double scale = 1., double scale_ico = 1.);
269 
274  inline const dealii::SmartPointer<const dealii::FESystem<dim> >
275  GetFESystem() const;
276 
281  inline const dealii::SmartPointer<
282  const dealii::hp::FECollection<dim> >
283  GetFECollection() const;
284 
289  inline std::string
290  GetDoFType() const;
291 
296  inline bool
297  HasFaces() const;
302  inline bool
303  HasPoints() const;
308  inline bool
309  HasInterfaces() const;
310 
315  inline dealii::UpdateFlags
316  GetUpdateFlags() const;
317 
322  inline dealii::UpdateFlags
323  GetFaceUpdateFlags() const;
324 
325  /******************************************************/
330  inline void
331  SetTime(double time,
332  unsigned int time_dof_number,
333  const TimeIterator& interval, bool initial = false);
334 
339  inline void
340  ComputeSparsityPattern(SPARSITYPATTERN & sparsity) const;
344  inline void
345  ComputeMGSparsityPattern(dealii::MGLevelObject<dealii::BlockSparsityPattern> & mg_sparsity_patterns,
346  unsigned int n_levels) const;
347 
348  inline void
349  ComputeMGSparsityPattern(dealii::MGLevelObject<dealii::SparsityPattern> & mg_sparsity_patterns,
350  unsigned int n_levels) const;
351 
356  inline const std::vector<unsigned int>&
357  GetDirichletColors() const;
362  inline const std::vector<bool>&
363  GetDirichletCompMask(unsigned int color) const;
368  inline const Function<dim>
369  &
370  GetDirichletValues(
371  unsigned int color,
372  const std::map<std::string, const dealii::Vector<double>*> &param_values,
373  const std::map<std::string, const VECTOR*> &domain_values) const;
378  inline const std::vector<unsigned int>&
379  GetBoundaryEquationColors() const;
384  inline const dealii::ConstraintMatrix&
385  GetDoFConstraints() const;
390  const dealii::Function<dim>&
391  GetInitialValues() const;
392  /******************************************************/
395  {
396  return opt_problem_.GetOutputHandler();
397  }
398  OPTPROBLEM&
400  {
401  return opt_problem_;
402  }
403 
404  template<typename ELEMENTITERATOR>
405  bool AtInterface(ELEMENTITERATOR& element, unsigned int face) const;
406 
407  protected:
408 
409  private:
410  PDE& pde_;
411  OPTPROBLEM& opt_problem_;
412 
413  std::vector<unsigned int> dirichlet_colors_;
414  std::vector<std::vector<bool> > dirichlet_comps_;
415  std::vector<PrimalDirichletData<DD, VECTOR, dim>*>
416  primal_dirichlet_values_;
417  std::vector<unsigned int> state_boundary_equation_colors_;
418  double interval_length_;
419  };
420 
421  /*****************************************************************************************/
422  /********************************IMPLEMENTATION*******************************************/
423  /*****************************************************************************************/
424 
425  template<typename OPTPROBLEM, typename PDE, typename DD,
426  typename SPARSITYPATTERN, typename VECTOR, int dim>
427  template<typename EDC>
428  void
429  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
430  dim>::ElementEquation(const EDC& edc,
431  dealii::Vector<double> &local_vector, double scale,
432  double scale_ico)
433  {
434  pde_.ElementEquation(edc, local_vector, scale*interval_length_, scale_ico*interval_length_);
435  }
436 
437  /******************************************************/
438 
439  template<typename OPTPROBLEM, typename PDE, typename DD,
440  typename SPARSITYPATTERN, typename VECTOR, int dim>
441  template<typename EDC>
442  void
443  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
444  dim>::ElementTimeEquation(const EDC& edc,
445  dealii::Vector<double> &local_vector, double scale)
446  {
447  pde_.ElementTimeEquation(edc, local_vector, scale);
448  }
449 
450  /******************************************************/
451 
452  template<typename OPTPROBLEM, typename PDE, typename DD,
453  typename SPARSITYPATTERN, typename VECTOR, int dim>
454  template<typename EDC>
455  void
456  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
457  dim>::ElementTimeEquationExplicit(const EDC& edc,
458  dealii::Vector<double> &local_vector, double scale)
459  {
460  pde_.ElementTimeEquationExplicit(edc, local_vector, scale);
461  }
462 
463  /******************************************************/
464 
465  template<typename OPTPROBLEM, typename PDE, typename DD,
466  typename SPARSITYPATTERN, typename VECTOR, int dim>
467  template<typename FDC>
468  void
469  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
470  dim>::FaceEquation(const FDC& fdc,
471  dealii::Vector<double> &local_vector, double scale, double scale_ico)
472  {
473  pde_.FaceEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
474  }
475 
476  /******************************************************/
477 
478  template<typename OPTPROBLEM, typename PDE, typename DD,
479  typename SPARSITYPATTERN, typename VECTOR, int dim>
480  template<typename FDC>
481  void
482  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
483  dim>::InterfaceEquation(const FDC& fdc,
484  dealii::Vector<double> &local_vector, double scale, double scale_ico)
485  {
486  pde_.InterfaceEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
487  }
488  /******************************************************/
489 
490  template<typename OPTPROBLEM, typename PDE, typename DD,
491  typename SPARSITYPATTERN, typename VECTOR, int dim>
492  template<typename FDC>
493  void
494  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
495  dim>::BoundaryEquation(const FDC& fdc,
496  dealii::Vector<double> &local_vector, double scale, double scale_ico)
497  {
498  pde_.BoundaryEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
499  }
500 
501  /******************************************************/
502 
503  template<typename OPTPROBLEM, typename PDE, typename DD,
504  typename SPARSITYPATTERN, typename VECTOR, int dim>
505  template<typename EDC>
506  void
507  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
508  dim>::ElementRhs(const EDC& edc,
509  dealii::Vector<double> &local_vector, double scale)
510  {
511  pde_.ElementRightHandSide(edc, local_vector, scale*interval_length_);
512  }
513 
514  /******************************************************/
515 
516  template<typename OPTPROBLEM, typename PDE, typename DD,
517  typename SPARSITYPATTERN, typename VECTOR, int dim>
518  void
520  const std::map<std::string, const dealii::Vector<double>*> &/*param_values*/,
521  const std::map<std::string, const VECTOR*> &/*domain_values*/,
522  VECTOR& /*rhs_vector*/, double /*scale*/)
523  {
524 
525  }
526 
527  /******************************************************/
528 
529  template<typename OPTPROBLEM, typename PDE, typename DD,
530  typename SPARSITYPATTERN, typename VECTOR, int dim>
531  template<typename FDC>
532  void
533  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
534  dim>::FaceRhs(const FDC& fdc,
535  dealii::Vector<double> &local_vector, double scale)
536  {
537  pde_.FaceRightHandSide(fdc, local_vector, scale*interval_length_);
538  }
539 
540  /******************************************************/
541 
542  template<typename OPTPROBLEM, typename PDE, typename DD,
543  typename SPARSITYPATTERN, typename VECTOR, int dim>
544  template<typename FDC>
545  void
546  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
547  dim>::BoundaryRhs(const FDC& fdc,
548  dealii::Vector<double> &local_vector, double scale)
549  {
550  pde_.BoundaryRightHandSide(fdc, local_vector, scale*interval_length_);
551  }
552 
553  /******************************************************/
554 
555  template<typename OPTPROBLEM, typename PDE, typename DD,
556  typename SPARSITYPATTERN, typename VECTOR, int dim>
557  template<typename EDC>
558  void
559  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
560  dim>::ElementMatrix(const EDC& edc,
561  dealii::FullMatrix<double> &local_entry_matrix, double scale,
562  double scale_ico)
563  {
564  pde_.ElementMatrix(edc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
565  }
566 
567  /******************************************************/
568 
569  template<typename OPTPROBLEM, typename PDE, typename DD,
570  typename SPARSITYPATTERN, typename VECTOR, int dim>
571  template<typename EDC>
572  void
573  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
574  dim>::ElementTimeMatrix(const EDC& edc,
575  FullMatrix<double> &local_entry_matrix)
576  {
577  pde_.ElementTimeMatrix(edc, local_entry_matrix);
578  }
579 
580  /******************************************************/
581 
582  template<typename OPTPROBLEM, typename PDE, typename DD,
583  typename SPARSITYPATTERN, typename VECTOR, int dim>
584  template<typename EDC>
585  void
586  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
587  dim>::ElementTimeMatrixExplicit(const EDC& edc,
588  dealii::FullMatrix<double> &local_entry_matrix)
589  {
590  pde_.ElementTimeMatrixExplicit(edc, local_entry_matrix);
591  }
592 
593  /******************************************************/
594 
595  template<typename OPTPROBLEM, typename PDE, typename DD,
596  typename SPARSITYPATTERN, typename VECTOR, int dim>
597  template<typename FDC>
598  void
599  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
600  dim>::FaceMatrix(const FDC& fdc,
601  FullMatrix<double> &local_entry_matrix, double scale,
602  double scale_ico)
603  {
604  pde_.FaceMatrix(fdc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
605  }
606 
607  /******************************************************/
608 
609  template<typename OPTPROBLEM, typename PDE, typename DD,
610  typename SPARSITYPATTERN, typename VECTOR, int dim>
611  template<typename FDC>
612  void
613  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
614  dim>::InterfaceMatrix(const FDC& fdc,
615  FullMatrix<double> &local_entry_matrix, double scale,
616  double scale_ico)
617  {
618  pde_.InterfaceMatrix(fdc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
619  }
620 
621  /******************************************************/
622 
623  template<typename OPTPROBLEM, typename PDE, typename DD,
624  typename SPARSITYPATTERN, typename VECTOR, int dim>
625  template<typename FDC>
626  void
627  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
628  dim>::BoundaryMatrix(const FDC& fdc,
629  FullMatrix<double> &local_matrix, double scale,
630  double scale_ico)
631  {
632  pde_.BoundaryMatrix(fdc, local_matrix, scale*interval_length_, scale_ico*interval_length_);
633  }
634 
635  /******************************************************/
636 
637  template<typename OPTPROBLEM, typename PDE, typename DD,
638  typename SPARSITYPATTERN, typename VECTOR, int dim>
639  std::string
641  {
642  return "state";
643  }
644 
645  /******************************************************/
646 
647  template<typename OPTPROBLEM, typename PDE, typename DD,
648  typename SPARSITYPATTERN, typename VECTOR, int dim>
649  const SmartPointer<const dealii::FESystem<dim> >
651  {
652  return opt_problem_.GetSpaceTimeHandler()->GetFESystem("state");
653  }
654 
655  /******************************************************/
656  template<typename OPTPROBLEM, typename PDE, typename DD,
657  typename SPARSITYPATTERN, typename VECTOR, int dim>
658  const SmartPointer<const dealii::hp::FECollection<dim> >
660  {
661  return opt_problem_.GetSpaceTimeHandler()->GetFECollection("state");
662  }
663 
664  /******************************************************/
665 
666  template<typename OPTPROBLEM, typename PDE, typename DD,
667  typename SPARSITYPATTERN, typename VECTOR, int dim>
668  UpdateFlags
670  {
671  UpdateFlags r;
672  r = pde_.GetUpdateFlags();
673  return r | update_JxW_values;
674  }
675 
676  /******************************************************/
677 
678  template<typename OPTPROBLEM, typename PDE, typename DD,
679  typename SPARSITYPATTERN, typename VECTOR, int dim>
680  UpdateFlags
682  {
683  UpdateFlags r;
684  r = pde_.GetFaceUpdateFlags();
685  return r | update_JxW_values;
686  }
687 
688  /******************************************************/
689 
690  template<typename OPTPROBLEM, typename PDE, typename DD,
691  typename SPARSITYPATTERN, typename VECTOR, int dim>
692  void
694  double time,
695  unsigned int time_dof_number, const TimeIterator& interval, bool initial)
696  {
697  opt_problem_.SetTime(time, time_dof_number, interval, initial);
698  interval_length_ = opt_problem_.GetSpaceTimeHandler()->GetStepSize();
699  }
700 
701  /******************************************************/
702 
703  template<typename OPTPROBLEM, typename PDE, typename DD,
704  typename SPARSITYPATTERN, typename VECTOR, int dim>
705  void
707  SPARSITYPATTERN & sparsity) const
708  {
709  opt_problem_.GetSpaceTimeHandler()->ComputeStateSparsityPattern(sparsity);
710  }
711 
712  /******************************************************/
713 
714  template<typename OPTPROBLEM, typename PDE, typename DD,
715  typename SPARSITYPATTERN, typename VECTOR, int dim>
716  void
718  dealii::MGLevelObject<dealii::BlockSparsityPattern> & mg_sparsity_patterns,
719  unsigned int n_levels) const
720  {
721  opt_problem_.GetSpaceTimeHandler()->ComputeMGStateSparsityPattern(mg_sparsity_patterns, n_levels);
722  }
723 
724 /******************************************************/
725 
726  template<typename OPTPROBLEM, typename PDE, typename DD,
727  typename SPARSITYPATTERN, typename VECTOR, int dim>
728  void
730  dealii::MGLevelObject<dealii::SparsityPattern> & mg_sparsity_patterns,
731  unsigned int n_levels) const
732  {
733  opt_problem_.GetSpaceTimeHandler()->ComputeMGStateSparsityPattern(mg_sparsity_patterns, n_levels);
734  }
735 
736 
737 
738  /******************************************************/
739 
740  template<typename OPTPROBLEM, typename PDE, typename DD,
741  typename SPARSITYPATTERN, typename VECTOR, int dim>
742  bool
744  {
745  return pde_.HasFaces();
746  }
747 
748  /******************************************************/
749 
750  template<typename OPTPROBLEM, typename PDE, typename DD,
751  typename SPARSITYPATTERN, typename VECTOR, int dim>
752  bool
754  {
755  return false;//We have no PointRhs in normal stateproblems at the moment.
756  }
757 
758 
759  /******************************************************/
760 
761  template<typename OPTPROBLEM, typename PDE, typename DD,
762  typename SPARSITYPATTERN, typename VECTOR, int dim>
763  bool
765  {
766  return pde_.HasInterfaces();
767  }
768 
769  /******************************************************/
770 
771  template<typename OPTPROBLEM, typename PDE, typename DD,
772  typename SPARSITYPATTERN, typename VECTOR, int dim>
773  const std::vector<unsigned int>&
775  {
776  return dirichlet_colors_;
777  }
778 
779  /******************************************************/
780 
781  template<typename OPTPROBLEM, typename PDE, typename DD,
782  typename SPARSITYPATTERN, typename VECTOR, int dim>
783  const std::vector<bool>&
785  unsigned int color) const
786  {
787  unsigned int comp = dirichlet_colors_.size();
788  for (unsigned int i = 0; i < dirichlet_colors_.size(); ++i)
789  {
790  if (dirichlet_colors_[i] == color)
791  {
792  comp = i;
793  break;
794  }
795  }
796  if (comp == dirichlet_colors_.size())
797  {
798  std::stringstream s;
799  s << "DirichletColor" << color << " has not been found !";
800  throw DOpEException(s.str(), "OptProblem::GetDirichletCompMask");
801  }
802  return dirichlet_comps_[comp];
803  }
804 
805  /******************************************************/
806 
807  template<typename OPTPROBLEM, typename PDE, typename DD,
808  typename SPARSITYPATTERN, typename VECTOR, int dim>
809  const Function<dim>&
811  unsigned int color,
812  const std::map<std::string, const dealii::Vector<double>*> &param_values,
813  const std::map<std::string, const VECTOR*> &domain_values) const
814  {
815  unsigned int col = dirichlet_colors_.size();
816  for (unsigned int i = 0; i < dirichlet_colors_.size(); ++i)
817  {
818  if (dirichlet_colors_[i] == color)
819  {
820  col = i;
821  break;
822  }
823  }
824  if (col == dirichlet_colors_.size())
825  {
826  std::stringstream s;
827  s << "DirichletColor" << color << " has not been found !";
828  throw DOpEException(s.str(), "OptProblem::GetDirichletValues");
829  }
830  primal_dirichlet_values_[col]->ReInit(param_values, domain_values, color);
831  return *(primal_dirichlet_values_[col]);
832  }
833 
834  /******************************************************/
835 
836  template<typename OPTPROBLEM, typename PDE, typename DD,
837  typename SPARSITYPATTERN, typename VECTOR, int dim>
838  const std::vector<unsigned int>&
840  {
841  return state_boundary_equation_colors_;
842  }
843 
844  /******************************************************/
845 
846  template<typename OPTPROBLEM, typename PDE, typename DD,
847  typename SPARSITYPATTERN, typename VECTOR, int dim>
848  const dealii::ConstraintMatrix&
850  {
851  return opt_problem_.GetSpaceTimeHandler()->GetStateDoFConstraints();
852  }
853 
854  /******************************************************/
855 
856  template<typename OPTPROBLEM, typename PDE, typename DD,
857  typename SPARSITYPATTERN, typename VECTOR, int dim> const dealii::Function<dim>&
859  {
860  return opt_problem_.GetInitialValues();
861  }
862 
863  /******************************************************/
864 
865  template<typename OPTPROBLEM, typename PDE, typename DD,
866  typename SPARSITYPATTERN, typename VECTOR, int dim>
867  template<typename ELEMENTITERATOR>
869  ::AtInterface(ELEMENTITERATOR& element, unsigned int face) const
870  {
871  return pde_.AtInterface(element,face);
872  }
873 
874 
876 }
877 #endif
Definition: stateproblem.h:48
void FaceMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_entry_matrix, double scale=1., double scale_ico=1.)
void InterfaceEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:483
dealii::UpdateFlags GetFaceUpdateFlags() const
Definition: stateproblem.h:681
void Init_ElementEquation(const EDC &edc, dealii::Vector< double > &local_vector, double scale, double scale_ico)
Definition: stateproblem.h:80
OPTPROBLEM & GetBaseProblem()
Definition: stateproblem.h:399
void ElementTimeEquationExplicit(const EDC &edc, dealii::Vector< double > &local_vector, double scale=1.)
Definition: stateproblem.h:457
Definition: dopetypes.h:106
void BoundaryEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:495
StateProblem(OPTPROBLEM &OP, PDE &pde)
Definition: stateproblem.h:51
const dealii::SmartPointer< const dealii::FESystem< dim > > GetFESystem() const
Definition: stateproblem.h:650
void SetTime(double time, unsigned int time_dof_number, const TimeIterator &interval, bool initial=false)
Definition: stateproblem.h:693
Definition: optproblemcontainer.h:70
void Init_PointRhs(const std::map< std::string, const dealii::Vector< double > * > &, const std::map< std::string, const VECTOR * > &, VECTOR &, double)
Definition: stateproblem.h:104
std::string GetName() const
Definition: stateproblem.h:63
const dealii::SmartPointer< const dealii::hp::FECollection< dim > > GetFECollection() const
Definition: stateproblem.h:659
void FaceEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:470
DOpEOutputHandler< VECTOR > * GetOutputHandler()
Definition: stateproblem.h:394
Definition: timeiterator.h:62
void ElementMatrix(const EDC &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:560
bool HasInterfaces() const
Definition: stateproblem.h:764
void ElementTimeEquation(const EDC &edc, dealii::Vector< double > &local_vector, double scale=1.)
Definition: stateproblem.h:444
dealii::UpdateFlags GetUpdateFlags() const
Definition: stateproblem.h:669
void Init_ElementRhs(const EDC &edc, dealii::Vector< double > &local_vector, double scale)
Definition: stateproblem.h:93
bool HasFaces() const
Definition: stateproblem.h:743
void InterfaceMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_entry_matrix, double scale=1., double scale_ico=1.)
void BoundaryMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_matrix, double scale=1., double scale_ico=1.)
std::string GetType() const
Definition: stateproblem.h:68
void ElementTimeMatrixExplicit(const EDC &edc, dealii::FullMatrix< double > &local_entry_matrix)
Definition: stateproblem.h:587
void Init_ElementMatrix(const EDC &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale, double scale_ico)
Definition: stateproblem.h:118
void ElementTimeMatrix(const EDC &edc, dealii::FullMatrix< double > &local_entry_matrix)
const dealii::Function< dim > & GetInitialValues() const
Definition: stateproblem.h:858
Definition: dopeexception.h:35
bool AtInterface(ELEMENTITERATOR &element, unsigned int face) const
Definition: stateproblem.h:869