Reference documentation for deal.II version 8.5.1
theta_timestepping.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii__theta_timestepping_h
18 #define dealii__theta_timestepping_h
19 
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/algorithms/operator.h>
22 #include <deal.II/algorithms/timestep_control.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 class ParameterHandler;
27 
28 namespace Algorithms
29 {
39  struct TimestepData
40  {
42  double time;
44  double step;
45  };
46 
189  template <typename VectorType>
191  {
192  public:
200 
212  virtual void operator() (AnyData &out, const AnyData &in);
213 
217  virtual void notify(const Event &);
218 
224 
228  static void declare_parameters (ParameterHandler &param);
229 
233  void parse_parameters (ParameterHandler &param);
234 
238  double current_time() const;
242  double step_size() const;
246  double theta() const;
247 
251  double theta(double new_theta);
252 
259  const TimestepData &explicit_data() const;
260 
267  const TimestepData &implicit_data() const;
268 
273 
274  private:
280 
285  double vtheta;
289  bool adaptive;
290 
295 
300 
301 
313 
325 
330  };
331 
332 
333  template <typename VectorType>
334  inline
335  const TimestepData &
337  {
338  return d_explicit;
339  }
340 
341 
342  template <typename VectorType>
343  inline
344  const TimestepData &
346  {
347  return d_implicit;
348  }
349 
350 
351  template <typename VectorType>
352  inline
355  {
356  return control;
357  }
358 
359  template <typename VectorType>
360  inline
362  {
363  output = &out;
364  }
365 
366 
367  template <typename VectorType>
368  inline
370  {
371  return vtheta;
372  }
373 
374 
375  template <typename VectorType>
376  inline
377  double ThetaTimestepping<VectorType>::theta (double new_theta)
378  {
379  const double tmp = vtheta;
380  vtheta = new_theta;
381  return tmp;
382  }
383 
384 
385  template <typename VectorType>
386  inline
388  {
389  return control.now();
390  }
391 }
392 
393 DEAL_II_NAMESPACE_CLOSE
394 
395 #endif
virtual void notify(const Event &)
ThetaTimestepping(OperatorBase &op_explicit, OperatorBase &op_implicit)
TimestepControl & timestep_control()
static void declare_parameters(ParameterHandler &param)
virtual void operator()(AnyData &out, const AnyData &in)
const TimestepData & explicit_data() const
void parse_parameters(ParameterHandler &param)
double time
The current time.
SmartPointer< OperatorBase, ThetaTimestepping< VectorType > > op_implicit
SmartPointer< OutputOperator< VectorType >, ThetaTimestepping< VectorType > > output
void set_output(OutputOperator< VectorType > &output)
double step
The current step size times something.
SmartPointer< OperatorBase, ThetaTimestepping< VectorType > > op_explicit
const TimestepData & implicit_data() const