www.mooseframework.org
AB2PredictorCorrector.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 // MOOSE includes
13 #include "TimeStepper.h"
14 
15 // C++ includes
16 #include <fstream>
17 
18 namespace libMesh
19 {
20 template <typename T>
21 class NumericVector;
22 }
23 
30 {
31 public:
33 
35 
36  virtual void step() override;
37  virtual void preExecute() override;
38  virtual void preSolve() override;
39  virtual void postSolve() override;
40  virtual bool converged() const override;
41 
42 protected:
43  virtual Real computeDT() override;
44  virtual Real computeInitialDT() override;
45 
47 
51 
54 
73  std::ofstream myfile;
74 };
AB2PredictorCorrector(const InputParameters &parameters)
Real & _error
global relative time discretization error estimate
virtual Real computeInitialDT() override
Called to compute _current_dt for the first timestep.
int _steps_between_increase
steps to take before increasing dt
Real _e_tol
error tolerance
NumericVector< Number > & _aux1
Base class for time stepping.
Definition: TimeStepper.h:22
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
NumericVector< Number > & _pred1
virtual void preSolve() override
Real & _dt_full
dt of the big step
virtual void step() override
Take a time step.
virtual void postSolve() override
virtual Real estimateTimeError(NumericVector< Number > &sol)
Real _e_max
maximal error
Real _max_increase
maximum increase ratio
Real _scaling_parameter
scaling_parameter for time step selection, default is 0.8
NumericVector< Number > & _u1
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool converged() const override
If the time step converged.
Real & _infnorm
infinity norm of the solution vector
const InputParameters & parameters() const
Get the parameters of the object.
virtual Real computeDT() override
Called to compute _current_dt for a normal step.
virtual void preExecute() override
int & _dt_steps_taken
steps taken at current dt
A TimeStepper based on the AB2 method.
static InputParameters validParams()