DOpE
backward_euler_problem.h
Go to the documentation of this file.
1 
24 #ifndef _BackwardEulerProblem_H_
25 #define _BackwardEulerProblem_H_
26 
27 #include "initialproblem.h"
28 #include "primal_ts_base.h"
29 
30 namespace DOpE
31 {
50  template<typename OPTPROBLEM, typename SPARSITYPATTERN, typename VECTOR,
51  int dopedim, int dealdim, template<int, int> class FE = dealii::FESystem,
52  template<int, int> class DH = dealii::DoFHandler>
53  class BackwardEulerProblem : public PrimalTSBase<OPTPROBLEM,
54  SPARSITYPATTERN, VECTOR, dopedim, dealdim, FE, DH>
55  {
56  public:
57  BackwardEulerProblem(OPTPROBLEM& OP) :
58  PrimalTSBase<OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim, dealdim,
59  FE, DH>(OP)
60  {
61  _initial_problem = NULL;
62  }
64  {
65  if (_initial_problem != NULL)
66  delete _initial_problem;
67  }
68 
69  /******************************************************/
70 
76  std::string
78  {
79  return "backward Euler";
80  }
81  /******************************************************/
82 
88  BackwardEulerProblem<OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim,
89  dealdim, FE, DH>, VECTOR, dealdim>&
91  {
92  if (_initial_problem == NULL)
93  {
94  _initial_problem = new InitialProblem<
95  BackwardEulerProblem<OPTPROBLEM, SPARSITYPATTERN, VECTOR,
96  dopedim, dealdim, FE, DH>, VECTOR, dealdim>(*this);
97  }
98  return *_initial_problem;
99  }
100 
101  /******************************************************/
102 
109  BackwardEulerProblem<OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim,
110  dealdim, FE, DH> &
112  {
113  return *this;
114  }
115 
116  /******************************************************/
117 
149  template<typename EDC>
150  void
151  ElementEquation(const EDC & edc,
152  dealii::Vector<double> &local_vector, double scale, double /*scale_ico*/)
153  {
154  if (this->GetPart() == "New")
155  {
156  dealii::Vector<double> tmp(local_vector);
157  tmp = 0.0;
158  this->GetProblem().ElementEquation(edc, tmp,
159  scale
160  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
161  scale
162  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
163  local_vector += tmp;
164 
165  tmp = 0.0;
166  this->GetProblem().ElementTimeEquation(edc, tmp, scale);
167  local_vector += tmp;
168 
169  this->GetProblem().ElementTimeEquationExplicit(edc, local_vector,
170  scale);
171  }
172  else if (this->GetPart() == "Old")
173  {
174  this->GetProblem().ElementTimeEquation(edc, local_vector,
175  (-1) * scale);
176  }
177  else
178  {
179  abort();
180  }
181  }
182 
183  /******************************************************/
184 
204  template<typename EDC>
205  void
206  ElementRhs(const EDC & edc,
207  dealii::Vector<double> &local_vector, double scale)
208  {
209  if (this->GetPart() == "New")
210  {
211  this->GetProblem().ElementRhs(edc, local_vector,
212  scale
213  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
214  }
215  else if (this->GetPart() == "Old")
216  {
217  }
218  else
219  {
220  abort();
221  }
222  }
223 
224  /******************************************************/
225 
250  void
252  const std::map<std::string, const dealii::Vector<double>*> &param_values,
253  const std::map<std::string, const VECTOR*> &domain_values,
254  VECTOR& rhs_vector, double scale)
255  {
256  if (this->GetPart() == "New")
257  {
258  this->GetProblem().PointRhs(param_values, domain_values, rhs_vector,
259  scale
260  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
261  }
262  else if (this->GetPart() == "Old")
263  {
264  }
265  else
266  {
267  abort();
268  }
269  }
270 
271  /******************************************************/
272 
304  template<typename EDC>
305  void
306  ElementMatrix(const EDC & edc,
307  dealii::FullMatrix<double> &local_matrix)
308  {
309  assert(this->GetPart() == "New");
310  dealii::FullMatrix<double> m(local_matrix);
311 
312  this->GetProblem().ElementMatrix(edc, local_matrix,
313  1.
314  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
315  1.
316  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
317 
318  m = 0.;
319  this->GetProblem().ElementTimeMatrix(edc, m);
320  local_matrix.add(1.0, m);
321 
322  m = 0.;
323  this->GetProblem().ElementTimeMatrixExplicit(edc, m);
324  local_matrix.add(1.0, m);
325 
326  }
327 
328  /******************************************************/
329 
345  template<typename FDC>
346  void
347  FaceEquation(const FDC & fdc,
348  dealii::Vector<double> &local_vector, double scale,
349  double /*scale_ico*/)
350  {
351  if (this->GetPart() == "New")
352  {
353  this->GetProblem().FaceEquation(fdc, local_vector,
354  scale
355  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
356  scale
357  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
358  }
359  else if (this->GetPart() == "Old")
360  {
361  }
362  else
363  {
364  abort();
365  }
366  }
367 
368  /******************************************************/
369 
386  template<typename FDC>
387  void
388  InterfaceEquation(const FDC& fdc,
389  dealii::Vector<double> &local_vector, double scale,
390  double /*scale_ico*/)
391  {
392  if (this->GetPart() == "New")
393  {
394  this->GetProblem().InterfaceEquation(fdc, local_vector,
395  scale
396  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
397  scale
398  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
399  }
400  else if (this->GetPart() == "Old")
401  {
402  }
403  else
404  {
405  abort();
406  }
407  }
408 
409  /******************************************************/
410 
425  template<typename FDC>
426  void
427  FaceRhs(const FDC & fdc,
428  dealii::Vector<double> &local_vector, double scale = 1.)
429  {
430  this->GetProblem().FaceRhs(fdc, local_vector,
431  scale
432  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
433  }
434 
435  /******************************************************/
436 
448  template<typename FDC>
449  void
450  FaceMatrix(const FDC & fdc,
451  dealii::FullMatrix<double> &local_matrix)
452  {
453  assert(this->GetPart() == "New");
454  this->GetProblem().FaceMatrix(fdc, local_matrix,
455  1.
456  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
457  1.
458  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
459 
460  }
461 
462  /******************************************************/
463 
475  template<typename FDC>
476  void
477  InterfaceMatrix(const FDC& fdc,
478  dealii::FullMatrix<double> &local_matrix)
479  {
480  assert(this->GetPart() == "New");
481  this->GetProblem().InterfaceMatrix(fdc, local_matrix,
482  1.
483  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
484  1.
485  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
486 
487  }
488 
489  /******************************************************/
490 
506  template<typename FDC>
507  void
508  BoundaryEquation(const FDC & fdc,
509  dealii::Vector<double> &local_vector, double scale,
510  double /*scale_ico*/)
511  {
512  if (this->GetPart() == "New")
513  {
514  this->GetProblem().BoundaryEquation(fdc, local_vector,
515  scale
516  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
517  scale
518  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
519  }
520  else if (this->GetPart() == "Old")
521  {
522  }
523  else
524  {
525  abort();
526  }
527 
528  }
529 
530  /******************************************************/
531 
545  template<typename FDC>
546  void
547  BoundaryRhs(const FDC & fdc,
548  dealii::Vector<double> &local_vector, double scale)
549  {
550  this->GetProblem().BoundaryRhs(fdc, local_vector,
551  scale
552  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
553  }
554 
555  /******************************************************/
556 
569  template<typename FDC>
570  void
571  BoundaryMatrix(const FDC & fdc,
572  dealii::FullMatrix<double> &local_matrix)
573  {
574  assert(this->GetPart() == "New");
575  this->GetProblem().BoundaryMatrix(fdc, local_matrix,
576  1.
577  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize(),
578  1.
579  * this->GetProblem().GetBaseProblem().GetSpaceTimeHandler()->GetStepSize());
580 
581  }
582  private:
584  BackwardEulerProblem<OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim,
585  dealdim, FE, DH>, VECTOR, dealdim> * _initial_problem;
586  };
587 }
588 
589 #endif
void FaceRhs(const FDC &fdc, dealii::Vector< double > &local_vector, double scale=1.)
Definition: backward_euler_problem.h:427
void ElementEquation(const EDC &edc, dealii::Vector< double > &local_vector, double scale, double)
Definition: backward_euler_problem.h:151
void InterfaceEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale, double)
Definition: backward_euler_problem.h:388
void PointRhs(const std::map< std::string, const dealii::Vector< double > * > &param_values, const std::map< std::string, const VECTOR * > &domain_values, VECTOR &rhs_vector, double scale)
Definition: backward_euler_problem.h:251
void BoundaryMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_matrix)
Definition: backward_euler_problem.h:571
void ElementRhs(const EDC &edc, dealii::Vector< double > &local_vector, double scale)
Definition: backward_euler_problem.h:206
BackwardEulerProblem(OPTPROBLEM &OP)
Definition: backward_euler_problem.h:57
Definition: primal_ts_base.h:49
Definition: backward_euler_problem.h:53
void BoundaryEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale, double)
Definition: backward_euler_problem.h:508
InitialProblem< BackwardEulerProblem< OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim, dealdim, FE, DH >, VECTOR, dealdim > & GetInitialProblem()
Definition: backward_euler_problem.h:90
BackwardEulerProblem< OPTPROBLEM, SPARSITYPATTERN, VECTOR, dopedim, dealdim, FE, DH > & GetBaseProblem()
Definition: backward_euler_problem.h:111
void BoundaryRhs(const FDC &fdc, dealii::Vector< double > &local_vector, double scale)
Definition: backward_euler_problem.h:547
~BackwardEulerProblem()
Definition: backward_euler_problem.h:63
Definition: initialproblem.h:42
const SpaceTimeHandler< FE, DH, SPARSITYPATTERN, VECTOR, dopedim, dealdim > * GetSpaceTimeHandler() const
Definition: ts_base.h:442
OPTPROBLEM & GetProblem()
Definition: ts_base.h:466
std::string GetName()
Definition: backward_euler_problem.h:77
std::string GetPart() const
Definition: ts_base.h:479
void FaceEquation(const FDC &fdc, dealii::Vector< double > &local_vector, double scale, double)
Definition: backward_euler_problem.h:347
void FaceMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_matrix)
Definition: backward_euler_problem.h:450
void InterfaceMatrix(const FDC &fdc, dealii::FullMatrix< double > &local_matrix)
Definition: backward_euler_problem.h:477
void ElementMatrix(const EDC &edc, dealii::FullMatrix< double > &local_matrix)
Definition: backward_euler_problem.h:306