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, const TimeIterator& interval, bool initial = false);
332 
337  inline void
338  ComputeSparsityPattern(SPARSITYPATTERN & sparsity) const;
342  inline void
343  ComputeMGSparsityPattern(dealii::MGLevelObject<dealii::BlockSparsityPattern> & mg_sparsity_patterns,
344  unsigned int n_levels) const;
345 
346  inline void
347  ComputeMGSparsityPattern(dealii::MGLevelObject<dealii::SparsityPattern> & mg_sparsity_patterns,
348  unsigned int n_levels) const;
349 
354  inline const std::vector<unsigned int>&
355  GetDirichletColors() const;
360  inline const std::vector<bool>&
361  GetDirichletCompMask(unsigned int color) const;
366  inline const Function<dim>
367  &
368  GetDirichletValues(
369  unsigned int color,
370  const std::map<std::string, const dealii::Vector<double>*> &param_values,
371  const std::map<std::string, const VECTOR*> &domain_values) const;
376  inline const std::vector<unsigned int>&
377  GetBoundaryEquationColors() const;
382  inline const dealii::ConstraintMatrix&
383  GetDoFConstraints() const;
388  const dealii::Function<dim>&
389  GetInitialValues() const;
390  /******************************************************/
393  {
394  return opt_problem_.GetOutputHandler();
395  }
396  OPTPROBLEM&
398  {
399  return opt_problem_;
400  }
401  protected:
402 
403  private:
404  PDE& pde_;
405  OPTPROBLEM& opt_problem_;
406 
407  std::vector<unsigned int> dirichlet_colors_;
408  std::vector<std::vector<bool> > dirichlet_comps_;
409  std::vector<PrimalDirichletData<DD, VECTOR, dim>*>
410  primal_dirichlet_values_;
411  std::vector<unsigned int> state_boundary_equation_colors_;
412  double interval_length_;
413  };
414 
415  /*****************************************************************************************/
416  /********************************IMPLEMENTATION*******************************************/
417  /*****************************************************************************************/
418 
419  template<typename OPTPROBLEM, typename PDE, typename DD,
420  typename SPARSITYPATTERN, typename VECTOR, int dim>
421  template<typename EDC>
422  void
423  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
424  dim>::ElementEquation(const EDC& edc,
425  dealii::Vector<double> &local_vector, double scale,
426  double scale_ico)
427  {
428  pde_.ElementEquation(edc, local_vector, scale*interval_length_, scale_ico*interval_length_);
429  }
430 
431  /******************************************************/
432 
433  template<typename OPTPROBLEM, typename PDE, typename DD,
434  typename SPARSITYPATTERN, typename VECTOR, int dim>
435  template<typename EDC>
436  void
437  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
438  dim>::ElementTimeEquation(const EDC& edc,
439  dealii::Vector<double> &local_vector, double scale)
440  {
441  pde_.ElementTimeEquation(edc, local_vector, scale);
442  }
443 
444  /******************************************************/
445 
446  template<typename OPTPROBLEM, typename PDE, typename DD,
447  typename SPARSITYPATTERN, typename VECTOR, int dim>
448  template<typename EDC>
449  void
450  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
451  dim>::ElementTimeEquationExplicit(const EDC& edc,
452  dealii::Vector<double> &local_vector, double scale)
453  {
454  pde_.ElementTimeEquationExplicit(edc, local_vector, scale);
455  }
456 
457  /******************************************************/
458 
459  template<typename OPTPROBLEM, typename PDE, typename DD,
460  typename SPARSITYPATTERN, typename VECTOR, int dim>
461  template<typename FDC>
462  void
463  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
464  dim>::FaceEquation(const FDC& fdc,
465  dealii::Vector<double> &local_vector, double scale, double scale_ico)
466  {
467  pde_.FaceEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
468  }
469 
470  /******************************************************/
471 
472  template<typename OPTPROBLEM, typename PDE, typename DD,
473  typename SPARSITYPATTERN, typename VECTOR, int dim>
474  template<typename FDC>
475  void
476  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
477  dim>::InterfaceEquation(const FDC& fdc,
478  dealii::Vector<double> &local_vector, double scale, double scale_ico)
479  {
480  pde_.InterfaceEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
481  }
482  /******************************************************/
483 
484  template<typename OPTPROBLEM, typename PDE, typename DD,
485  typename SPARSITYPATTERN, typename VECTOR, int dim>
486  template<typename FDC>
487  void
488  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
489  dim>::BoundaryEquation(const FDC& fdc,
490  dealii::Vector<double> &local_vector, double scale, double scale_ico)
491  {
492  pde_.BoundaryEquation(fdc, local_vector, scale*interval_length_, scale_ico*interval_length_);
493  }
494 
495  /******************************************************/
496 
497  template<typename OPTPROBLEM, typename PDE, typename DD,
498  typename SPARSITYPATTERN, typename VECTOR, int dim>
499  template<typename EDC>
500  void
501  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
502  dim>::ElementRhs(const EDC& edc,
503  dealii::Vector<double> &local_vector, double scale)
504  {
505  pde_.ElementRightHandSide(edc, local_vector, scale*interval_length_);
506  }
507 
508  /******************************************************/
509 
510  template<typename OPTPROBLEM, typename PDE, typename DD,
511  typename SPARSITYPATTERN, typename VECTOR, int dim>
512  void
514  const std::map<std::string, const dealii::Vector<double>*> &/*param_values*/,
515  const std::map<std::string, const VECTOR*> &/*domain_values*/,
516  VECTOR& /*rhs_vector*/, double /*scale*/)
517  {
518 
519  }
520 
521  /******************************************************/
522 
523  template<typename OPTPROBLEM, typename PDE, typename DD,
524  typename SPARSITYPATTERN, typename VECTOR, int dim>
525  template<typename FDC>
526  void
527  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
528  dim>::FaceRhs(const FDC& fdc,
529  dealii::Vector<double> &local_vector, double scale)
530  {
531  pde_.FaceRightHandSide(fdc, local_vector, scale*interval_length_);
532  }
533 
534  /******************************************************/
535 
536  template<typename OPTPROBLEM, typename PDE, typename DD,
537  typename SPARSITYPATTERN, typename VECTOR, int dim>
538  template<typename FDC>
539  void
540  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
541  dim>::BoundaryRhs(const FDC& fdc,
542  dealii::Vector<double> &local_vector, double scale)
543  {
544  pde_.BoundaryRightHandSide(fdc, local_vector, scale*interval_length_);
545  }
546 
547  /******************************************************/
548 
549  template<typename OPTPROBLEM, typename PDE, typename DD,
550  typename SPARSITYPATTERN, typename VECTOR, int dim>
551  template<typename EDC>
552  void
553  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
554  dim>::ElementMatrix(const EDC& edc,
555  dealii::FullMatrix<double> &local_entry_matrix, double scale,
556  double scale_ico)
557  {
558  pde_.ElementMatrix(edc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
559  }
560 
561  /******************************************************/
562 
563  template<typename OPTPROBLEM, typename PDE, typename DD,
564  typename SPARSITYPATTERN, typename VECTOR, int dim>
565  template<typename EDC>
566  void
567  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
568  dim>::ElementTimeMatrix(const EDC& edc,
569  FullMatrix<double> &local_entry_matrix)
570  {
571  pde_.ElementTimeMatrix(edc, local_entry_matrix);
572  }
573 
574  /******************************************************/
575 
576  template<typename OPTPROBLEM, typename PDE, typename DD,
577  typename SPARSITYPATTERN, typename VECTOR, int dim>
578  template<typename EDC>
579  void
580  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
581  dim>::ElementTimeMatrixExplicit(const EDC& edc,
582  dealii::FullMatrix<double> &local_entry_matrix)
583  {
584  pde_.ElementTimeMatrixExplicit(edc, local_entry_matrix);
585  }
586 
587  /******************************************************/
588 
589  template<typename OPTPROBLEM, typename PDE, typename DD,
590  typename SPARSITYPATTERN, typename VECTOR, int dim>
591  template<typename FDC>
592  void
593  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
594  dim>::FaceMatrix(const FDC& fdc,
595  FullMatrix<double> &local_entry_matrix, double scale,
596  double scale_ico)
597  {
598  pde_.FaceMatrix(fdc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
599  }
600 
601  /******************************************************/
602 
603  template<typename OPTPROBLEM, typename PDE, typename DD,
604  typename SPARSITYPATTERN, typename VECTOR, int dim>
605  template<typename FDC>
606  void
607  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
608  dim>::InterfaceMatrix(const FDC& fdc,
609  FullMatrix<double> &local_entry_matrix, double scale,
610  double scale_ico)
611  {
612  pde_.InterfaceMatrix(fdc, local_entry_matrix, scale*interval_length_, scale_ico*interval_length_);
613  }
614 
615  /******************************************************/
616 
617  template<typename OPTPROBLEM, typename PDE, typename DD,
618  typename SPARSITYPATTERN, typename VECTOR, int dim>
619  template<typename FDC>
620  void
621  StateProblem<OPTPROBLEM, PDE, DD, SPARSITYPATTERN, VECTOR,
622  dim>::BoundaryMatrix(const FDC& fdc,
623  FullMatrix<double> &local_matrix, double scale,
624  double scale_ico)
625  {
626  pde_.BoundaryMatrix(fdc, local_matrix, scale*interval_length_, scale_ico*interval_length_);
627  }
628 
629  /******************************************************/
630 
631  template<typename OPTPROBLEM, typename PDE, typename DD,
632  typename SPARSITYPATTERN, typename VECTOR, int dim>
633  std::string
635  {
636  return "state";
637  }
638 
639  /******************************************************/
640 
641  template<typename OPTPROBLEM, typename PDE, typename DD,
642  typename SPARSITYPATTERN, typename VECTOR, int dim>
643  const SmartPointer<const dealii::FESystem<dim> >
645  {
646  return opt_problem_.GetSpaceTimeHandler()->GetFESystem("state");
647  }
648 
649  /******************************************************/
650  template<typename OPTPROBLEM, typename PDE, typename DD,
651  typename SPARSITYPATTERN, typename VECTOR, int dim>
652  const SmartPointer<const dealii::hp::FECollection<dim> >
654  {
655  return opt_problem_.GetSpaceTimeHandler()->GetFECollection("state");
656  }
657 
658  /******************************************************/
659 
660  template<typename OPTPROBLEM, typename PDE, typename DD,
661  typename SPARSITYPATTERN, typename VECTOR, int dim>
662  UpdateFlags
664  {
665  UpdateFlags r;
666  r = pde_.GetUpdateFlags();
667  return r | update_JxW_values;
668  }
669 
670  /******************************************************/
671 
672  template<typename OPTPROBLEM, typename PDE, typename DD,
673  typename SPARSITYPATTERN, typename VECTOR, int dim>
674  UpdateFlags
676  {
677  UpdateFlags r;
678  r = pde_.GetFaceUpdateFlags();
679  return r | update_JxW_values;
680  }
681 
682  /******************************************************/
683 
684  template<typename OPTPROBLEM, typename PDE, typename DD,
685  typename SPARSITYPATTERN, typename VECTOR, int dim>
686  void
688  double time, const TimeIterator& interval, bool initial)
689  {
690  opt_problem_.SetTime(time, interval, initial);
691  interval_length_ = opt_problem_.GetSpaceTimeHandler()->GetStepSize();
692  }
693 
694  /******************************************************/
695 
696  template<typename OPTPROBLEM, typename PDE, typename DD,
697  typename SPARSITYPATTERN, typename VECTOR, int dim>
698  void
700  SPARSITYPATTERN & sparsity) const
701  {
702  opt_problem_.GetSpaceTimeHandler()->ComputeStateSparsityPattern(sparsity);
703  }
704 
705  /******************************************************/
706 
707  template<typename OPTPROBLEM, typename PDE, typename DD,
708  typename SPARSITYPATTERN, typename VECTOR, int dim>
709  void
711  dealii::MGLevelObject<dealii::BlockSparsityPattern> & mg_sparsity_patterns,
712  unsigned int n_levels) const
713  {
714  opt_problem_.GetSpaceTimeHandler()->ComputeMGStateSparsityPattern(mg_sparsity_patterns, n_levels);
715  }
716 
717 /******************************************************/
718 
719  template<typename OPTPROBLEM, typename PDE, typename DD,
720  typename SPARSITYPATTERN, typename VECTOR, int dim>
721  void
723  dealii::MGLevelObject<dealii::SparsityPattern> & mg_sparsity_patterns,
724  unsigned int n_levels) const
725  {
726  opt_problem_.GetSpaceTimeHandler()->ComputeMGStateSparsityPattern(mg_sparsity_patterns, n_levels);
727  }
728 
729 
730 
731  /******************************************************/
732 
733  template<typename OPTPROBLEM, typename PDE, typename DD,
734  typename SPARSITYPATTERN, typename VECTOR, int dim>
735  bool
737  {
738  return pde_.HasFaces();
739  }
740 
741  /******************************************************/
742 
743  template<typename OPTPROBLEM, typename PDE, typename DD,
744  typename SPARSITYPATTERN, typename VECTOR, int dim>
745  bool
747  {
748  return false;//We have no PointRhs in normal stateproblems at the moment.
749  }
750 
751 
752  /******************************************************/
753 
754  template<typename OPTPROBLEM, typename PDE, typename DD,
755  typename SPARSITYPATTERN, typename VECTOR, int dim>
756  bool
758  {
759  return pde_.HasInterfaces();
760  }
761 
762  /******************************************************/
763 
764  template<typename OPTPROBLEM, typename PDE, typename DD,
765  typename SPARSITYPATTERN, typename VECTOR, int dim>
766  const std::vector<unsigned int>&
768  {
769  return dirichlet_colors_;
770  }
771 
772  /******************************************************/
773 
774  template<typename OPTPROBLEM, typename PDE, typename DD,
775  typename SPARSITYPATTERN, typename VECTOR, int dim>
776  const std::vector<bool>&
778  unsigned int color) const
779  {
780  unsigned int comp = dirichlet_colors_.size();
781  for (unsigned int i = 0; i < dirichlet_colors_.size(); ++i)
782  {
783  if (dirichlet_colors_[i] == color)
784  {
785  comp = i;
786  break;
787  }
788  }
789  if (comp == dirichlet_colors_.size())
790  {
791  std::stringstream s;
792  s << "DirichletColor" << color << " has not been found !";
793  throw DOpEException(s.str(), "OptProblem::GetDirichletCompMask");
794  }
795  return dirichlet_comps_[comp];
796  }
797 
798  /******************************************************/
799 
800  template<typename OPTPROBLEM, typename PDE, typename DD,
801  typename SPARSITYPATTERN, typename VECTOR, int dim>
802  const Function<dim>&
804  unsigned int color,
805  const std::map<std::string, const dealii::Vector<double>*> &param_values,
806  const std::map<std::string, const VECTOR*> &domain_values) const
807  {
808  unsigned int col = dirichlet_colors_.size();
809  for (unsigned int i = 0; i < dirichlet_colors_.size(); ++i)
810  {
811  if (dirichlet_colors_[i] == color)
812  {
813  col = i;
814  break;
815  }
816  }
817  if (col == dirichlet_colors_.size())
818  {
819  std::stringstream s;
820  s << "DirichletColor" << color << " has not been found !";
821  throw DOpEException(s.str(), "OptProblem::GetDirichletValues");
822  }
823  primal_dirichlet_values_[col]->ReInit(param_values, domain_values, color);
824  return *(primal_dirichlet_values_[col]);
825  }
826 
827  /******************************************************/
828 
829  template<typename OPTPROBLEM, typename PDE, typename DD,
830  typename SPARSITYPATTERN, typename VECTOR, int dim>
831  const std::vector<unsigned int>&
833  {
834  return state_boundary_equation_colors_;
835  }
836 
837  /******************************************************/
838 
839  template<typename OPTPROBLEM, typename PDE, typename DD,
840  typename SPARSITYPATTERN, typename VECTOR, int dim>
841  const dealii::ConstraintMatrix&
843  {
844  return opt_problem_.GetSpaceTimeHandler()->GetStateDoFConstraints();
845  }
846  template<typename OPTPROBLEM, typename PDE, typename DD,
847  typename SPARSITYPATTERN, typename VECTOR, int dim> const dealii::Function<dim>&
849  {
850  return opt_problem_.GetInitialValues();
851  }
853 }
854 #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:477
dealii::UpdateFlags GetFaceUpdateFlags() const
Definition: stateproblem.h:675
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:397
void ElementTimeEquationExplicit(const EDC &edc, dealii::Vector< double > &local_vector, double scale=1.)
Definition: stateproblem.h:451
Definition: dopetypes.h:106
void BoundaryEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:489
StateProblem(OPTPROBLEM &OP, PDE &pde)
Definition: stateproblem.h:51
const dealii::SmartPointer< const dealii::FESystem< dim > > GetFESystem() const
Definition: stateproblem.h:644
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:653
void FaceEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:464
DOpEOutputHandler< VECTOR > * GetOutputHandler()
Definition: stateproblem.h:392
Definition: timeiterator.h:63
void ElementMatrix(const EDC &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale=1., double scale_ico=1.)
Definition: stateproblem.h:554
bool HasInterfaces() const
Definition: stateproblem.h:757
void ElementTimeEquation(const EDC &edc, dealii::Vector< double > &local_vector, double scale=1.)
Definition: stateproblem.h:438
void SetTime(double time, const TimeIterator &interval, bool initial=false)
Definition: stateproblem.h:687
dealii::UpdateFlags GetUpdateFlags() const
Definition: stateproblem.h:663
void Init_ElementRhs(const EDC &edc, dealii::Vector< double > &local_vector, double scale)
Definition: stateproblem.h:93
bool HasFaces() const
Definition: stateproblem.h:736
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:581
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:848
Definition: dopeexception.h:35