DOpE
initialproblem.h
Go to the documentation of this file.
1 
24 #ifndef INITIAL_PROBLEM_H_
25 #define INITIAL_PROBLEM_H_
26 
27 #include "spacetimehandler.h"
28 
29 using namespace dealii;
30 
31 namespace DOpE
32 {
41  template<typename PDE, typename VECTOR, int dim>
43  {
44  public:
45  InitialProblem(PDE& pde) :
46  pde_(pde)
47  {
48  }
49 
50  std::string
51  GetName() const
52  {
53  return "InitialProblem";
54  }
55  std::string
56  GetType() const
57  {
58  return "initial_state";
59  }
60 
61  /******************************************************/
66  template<typename EDC>
67  inline void
68  ElementEquation(const EDC& edc,
69  dealii::Vector<double> &local_vector, double scale,
70  double scale_ico);
71 
76  template<typename EDC>
77  inline void
78  ElementRhs(const EDC& edc,
79  dealii::Vector<double> &local_vector, double scale = 1.);
80 
85  inline void
86  PointRhs(
87  const std::map<std::string, const dealii::Vector<double>*> &param_values,
88  const std::map<std::string, const VECTOR*> &domain_values,
89  VECTOR& rhs_vector, double scale = 1);
90 
95  template<typename EDC>
96  inline void
97  ElementMatrix(const EDC& edc,
98  dealii::FullMatrix<double> &local_matrix, double scale = 1.,
99  double scale_ico = 1.);
100 
105  template<typename FDC>
106  inline void
107  FaceEquation(const FDC& fdc,
108  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
109 
114  template<typename FDC>
115  inline void
116  InterfaceEquation(const FDC& fdc,
117  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
118 
123  template<typename FDC>
124  inline void
125  FaceRhs(const FDC& fdc,
126  dealii::Vector<double> &local_vector, double scale = 1.);
127 
132  template<typename FDC>
133  inline void
134  FaceMatrix(const FDC& fdc,
135  dealii::FullMatrix<double> &local_matrix, double scale = 1., double scale_ico = 1.);
136 
141  template<typename FDC>
142  inline void
143  InterfaceMatrix(const FDC& fdc,
144  dealii::FullMatrix<double> &local_matrix, double scale = 1., double scale_ico = 1.);
145 
150  template<typename FDC>
151  inline void
152  BoundaryEquation(const FDC& fdc,
153  dealii::Vector<double> &local_vector, double scale = 1., double scale_ico = 1.);
154 
159  template<typename FDC>
160  inline void
161  BoundaryRhs(const FDC& fdc,
162  dealii::Vector<double> &local_vector, double scale = 1.);
163 
168  template<typename FDC>
169  inline void
170  BoundaryMatrix(const FDC& fdc,
171  dealii::FullMatrix<double> &local_matrix, double scale = 1., double scale_ico = 1.);
172 
177  inline const dealii::SmartPointer<const dealii::FESystem<dim> >
178  GetFESystem() const;
179 
184  inline const dealii::SmartPointer<const dealii::hp::FECollection<dim> >
185  GetFECollection() const;
186 
191  inline std::string
192  GetDoFType() const;
193 
198  inline bool
199  HasFaces() const;
204  inline bool
205  HasPoints() const;
210  inline bool
211  HasInterfaces() const;
212 
217  inline dealii::UpdateFlags
218  GetUpdateFlags() const;
219 
224  inline dealii::UpdateFlags
225  GetFaceUpdateFlags() const;
226 
227  /******************************************************/
232  inline void
233  SetTime(double time, const TimeIterator& interval);
234 
239  template<typename SPARSITYPATTERN>
240  inline void
241  ComputeSparsityPattern(SPARSITYPATTERN & sparsity) const;
242 
247  inline const std::vector<unsigned int>&
248  GetDirichletColors() const;
253  inline const std::vector<bool>&
254  GetDirichletCompMask(unsigned int color) const;
259  inline const Function<dim>&
260  GetDirichletValues(unsigned int color,
261  const std::map<std::string, const dealii::Vector<double>*> &param_values,
262  const std::map<std::string, const VECTOR*> &domain_values) const;
267  inline const std::vector<unsigned int>&
268  GetBoundaryEquationColors() const;
273  inline const dealii::ConstraintMatrix&
274  GetDoFConstraints() const;
279  const dealii::Function<dim>&
280  GetInitialValues() const;
281  /******************************************************/
284  {
285  return pde_.GetOutputHandler();
286  }
287  PDE&
289  {
290  return pde_;
291  }
292  template<typename ELEMENTITERATOR>
293  bool AtInterface(ELEMENTITERATOR& element, unsigned int face) const;
294 
295  protected:
296 
297  private:
298  PDE& pde_;
299  };
300 
301  /*****************************************************************************************/
302  /********************************IMPLEMENTATION*******************************************/
303  /*****************************************************************************************/
304 
305  template<typename PDE, typename VECTOR, int dim>
306  template<typename EDC>
307  void
309  dealii::Vector<double> &local_vector, double scale,
310  double scale_ico)
311  {
312  pde_.Init_ElementEquation(edc, local_vector, scale, scale_ico);
313  }
314 
315  /******************************************************/
316 
317  template<typename PDE, typename VECTOR, int dim>
318  template<typename FDC>
319  void
321  const FDC& fdc,
322  dealii::Vector<double> &local_vector, double scale, double scale_ico)
323  {
324  pde_.Init_FaceEquation(fdc, local_vector, scale, scale_ico);
325  }
326 
327  /******************************************************/
328 
329  template<typename PDE, typename VECTOR, int dim>
330  template<typename FDC>
331  void
333  const FDC& fdc,
334  dealii::Vector<double> &local_vector, double scale, double scale_ico)
335  {
336  pde_.Init_InterfaceEquation(fdc, local_vector, scale, scale_ico);
337  }
338  /******************************************************/
339 
340  template<typename PDE, typename VECTOR, int dim>
341  template<typename FDC>
342  void
344  const FDC& fdc,
345  dealii::Vector<double> &local_vector, double scale, double scale_ico)
346  {
347  pde_.Init_BoundaryEquation(fdc, local_vector, scale, scale_ico);
348  }
349 
350  /******************************************************/
351 
352  template<typename PDE, typename VECTOR, int dim>
353  template<typename EDC>
354  void
356  dealii::Vector<double> &local_vector, double scale)
357  {
358  pde_.Init_ElementRhs(edc, local_vector, scale);
359  }
360 
361  /******************************************************/
362  template<typename PDE, typename VECTOR, int dim>
363  void
365  const std::map<std::string, const dealii::Vector<double>*> &param_values,
366  const std::map<std::string, const VECTOR*> &domain_values,
367  VECTOR& rhs_vector, double scale)
368  {
369  pde_.Init_PointRhs(param_values, domain_values, rhs_vector, scale);
370  }
371 
372  /******************************************************/
373 
374  template<typename PDE, typename VECTOR, int dim>
375  template<typename FDC>
376  void
378  const FDC& /*fdc*/,
379  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
380  {
381  }
382 
383  /******************************************************/
384 
385  template<typename PDE, typename VECTOR, int dim>
386  template<typename FDC>
387  void
389  const FDC& /*fdc*/,
390  dealii::Vector<double> &/*local_vector*/, double /*scale*/)
391  {
392  }
393 
394  /******************************************************/
395 
396  template<typename PDE, typename VECTOR, int dim>
397  template<typename EDC>
398  void
400  dealii::FullMatrix<double> &local_matrix, double scale, double scale_ico)
401  {
402  pde_.Init_ElementMatrix(edc, local_matrix, scale, scale_ico);
403  }
404 
405  /******************************************************/
406 
407  template<typename PDE, typename VECTOR, int dim>
408  template<typename FDC>
409  void
411  FullMatrix<double> &local_matrix, double scale, double scale_ico)
412  {
413  pde_.Init_FaceMatrix(fdc, local_matrix, scale, scale_ico);
414  }
415 
416  /******************************************************/
417 
418  template<typename PDE, typename VECTOR, int dim>
419  template<typename FDC>
420  void
422  const FDC& fdc, FullMatrix<double> &local_matrix, double scale, double scale_ico)
423  {
424  pde_.Init_InterfaceMatrix(fdc, local_matrix, scale, scale_ico);
425  }
426 
427  /******************************************************/
428 
429  template<typename PDE, typename VECTOR, int dim>
430  template<typename FDC>
431  void
433  const FDC& fdc, FullMatrix<double> &local_matrix, double scale, double scale_ico)
434  {
435  pde_.Init_BoundaryMatrix(fdc, local_matrix, scale, scale_ico);
436  }
437 
438  /******************************************************/
439 
440  template<typename PDE, typename VECTOR, int dim>
441  std::string
443  {
444  return pde_.GetDoFType();
445  }
446 
447  /******************************************************/
448 
449  template<typename PDE, typename VECTOR, int dim>
450  const SmartPointer<const dealii::FESystem<dim> >
452  {
453  return pde_.GetFESystem();
454  }
455 
456  /******************************************************/
457  template<typename PDE, typename VECTOR, int dim>
458  const SmartPointer<const dealii::hp::FECollection<dim> >
460  {
461  return pde_.GetFECollection();
462  }
463 
464  /******************************************************/
465 
466  template<typename PDE, typename VECTOR, int dim>
467  UpdateFlags
469  {
470  UpdateFlags r;
471  r = pde_.GetUpdateFlags();
472  return r | update_JxW_values;
473  }
474 
475  /******************************************************/
476 
477  template<typename PDE, typename VECTOR, int dim>
478  UpdateFlags
480  {
481  UpdateFlags r;
482  r = pde_.GetFaceUpdateFlags();
483  return r | update_JxW_values;
484  }
485 
486  /******************************************************/
487 
488  template<typename PDE, typename VECTOR, int dim>
489  void
491  const TimeIterator& interval)
492  {
493  pde_.SetTime(time, interval);
494  }
495 
496  /******************************************************/
497 
498  template<typename PDE, typename VECTOR, int dim>
499  template<typename SPARSITYPATTERN>
500  void
502  SPARSITYPATTERN & sparsity) const
503  {
504  pde_.ComputeStateSparsityPattern(sparsity);
505  }
506 
507  /******************************************************/
508 
509  template<typename PDE, typename VECTOR, int dim>
510  bool
512  {
513  return pde_.HasFaces();
514  }
515 
516  /******************************************************/
517 
518  template<typename PDE, typename VECTOR, int dim>
519  bool
521  {
522  return pde_.HasPoints();
523  }
524 
525  /******************************************************/
526 
527  template<typename PDE, typename VECTOR, int dim>
528  bool
530  {
531  return pde_.HasInterfaces();
532  }
533 
534  /******************************************************/
535 
536  template<typename PDE, typename VECTOR, int dim>
537  const std::vector<unsigned int>&
539  {
540  return pde_.GetDirichletColors();
541  }
542 
543  /******************************************************/
544 
545  template<typename PDE, typename VECTOR, int dim>
546  const std::vector<bool>&
548  unsigned int color) const
549  {
550  return pde_.GetDirichletCompMask(color);
551  }
552 
553  /******************************************************/
554 
555  template<typename PDE, typename VECTOR, int dim>
556  const Function<dim>&
558  const std::map<std::string, const dealii::Vector<double>*> &param_values,
559  const std::map<std::string, const VECTOR*> &domain_values) const
560  {
561  return pde_.GetDirichletValues(color, param_values, domain_values);
562  }
563 
564  /******************************************************/
565 
566  template<typename PDE, typename VECTOR, int dim>
567  const std::vector<unsigned int>&
569  {
570  return pde_.GetBoundaryEquationColors();
571  }
572 
573  /******************************************************/
574 
575  template<typename PDE, typename VECTOR, int dim>
576  const dealii::ConstraintMatrix&
578  {
579  return pde_.GetDoFConstraints();
580  }
581 
582  /******************************************************/
583 
584  template<typename PDE, typename VECTOR, int dim>
585  const dealii::Function<dim>&
587  {
588  return pde_.GetInitialValues();
589  }
590 
591  /******************************************************/
592 
593  template<typename PDE, typename VECTOR, int dim>
594  template<typename ELEMENTITERATOR>
596  ::AtInterface(ELEMENTITERATOR& element, unsigned int face) const
597  {
598  return pde_.AtInterface(element,face);
599  }
600 
602 }
603 #endif
InitialProblem(PDE &pde)
Definition: initialproblem.h:45
dealii::UpdateFlags GetUpdateFlags() const
Definition: initialproblem.h:468
Definition: optproblemcontainer.h:70
std::string GetName() const
Definition: initialproblem.h:51
bool HasInterfaces() const
Definition: initialproblem.h:529
std::string GetType() const
Definition: initialproblem.h:56
std::string GetDoFType() const
Definition: initialproblem.h:442
Definition: timeiterator.h:62
bool HasPoints() const
Definition: initialproblem.h:520
const std::vector< unsigned int > & GetBoundaryEquationColors() const
Definition: initialproblem.h:568
Definition: initialproblem.h:42
bool AtInterface(ELEMENTITERATOR &element, unsigned int face) const
Definition: initialproblem.h:596
const std::vector< bool > & GetDirichletCompMask(unsigned int color) const
Definition: initialproblem.h:547
PDE & GetBaseProblem()
Definition: initialproblem.h:288
const dealii::Function< dim > & GetInitialValues() const
Definition: initialproblem.h:586
void SetTime(double time, const TimeIterator &interval)
Definition: initialproblem.h:490
bool HasFaces() const
Definition: initialproblem.h:511
const Function< dim > & GetDirichletValues(unsigned int color, const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values) const
Definition: initialproblem.h:557
const dealii::SmartPointer< const dealii::hp::FECollection< dim > > GetFECollection() const
Definition: initialproblem.h:459
const dealii::ConstraintMatrix & GetDoFConstraints() const
Definition: initialproblem.h:577
const dealii::SmartPointer< const dealii::FESystem< dim > > GetFESystem() const
Definition: initialproblem.h:451
DOpEOutputHandler< VECTOR > * GetOutputHandler()
Definition: initialproblem.h:283
dealii::UpdateFlags GetFaceUpdateFlags() const
Definition: initialproblem.h:479
const std::vector< unsigned int > & GetDirichletColors() const
Definition: initialproblem.h:538