timedisc_base_mod Module Reference

Data Types

interface  Finalize
 
interface  SolveODE
 
type  timedisc_base
 

Functions/Subroutines

subroutine setoutput (this, Mesh, Physics, config, IO)
 
subroutine acceptsolution (this, Mesh, Physics, Sources, Fluxes, time, dt, iter)
 
subroutine checkdata (this, Mesh, Physics, Fluxes, pvar, cvar, checkdatabm)
 compute the RHS of the spatially discretized PDE More...
 
pure integer function getorder (this)
 

Variables

integer, parameter, public modified_euler = 1
 
integer, parameter, public rk_fehlberg = 2
 
integer, parameter, public cash_karp = 3
 
integer, parameter, public dormand_prince = 4
 
integer, parameter, public ssprk = 5
 
integer, parameter, public dtcause_cfl = 0
 
integer, parameter, public dtcause_erradj = -1
 
integer, parameter, public dtcause_smallerr = -2
 
integer, parameter, public dtcause_fileio = -4
 smallest ts due to fileio More...
 
integer, parameter, public check_nothing = INT(b'0')
 
integer, parameter, public check_all = NOT(CHECK_NOTHING)
 
integer, parameter, public check_csound = INT(b'1')
 
integer, parameter, public check_pmin = INT(b'10')
 
integer, parameter, public check_rhomin = INT(b'100')
 
integer, parameter, public check_invalid = INT(b'1000')
 
integer, parameter, public check_tmin = INT(b'10000')
 
character(len=40), dimension(3), parameter fargo_method = (/ "dynamic background velocity ", "user supplied fixed background velocity ", "shearing box " /)
 
subroutine inittimedisc (this, Mesh, Physics, config, IO, ttype, tname)
 
subroutine integrationstep (this, Mesh, Physics, Sources, Fluxes, iter, config, IO)
 
real function adjusttimestep (this, maxerr, dtold)
 adjust the time step More...
 
real function calctimestep (this, Mesh, Physics, Sources, Fluxes, time, dtcause)
 Determines the CFL time step and time step limits due to source terms. More...
 
subroutine rejectsolution (this, Mesh, Physics, Sources, Fluxes, time, dt)
 
real function computeerror (this, Mesh, Physics, cvar_high, cvar_low)
 
subroutine computerhs (this, Mesh, Physics, Sources, Fluxes, time, dt, pvar, cvar, checkdatabm, rhs)
 compute the RHS of the spatially discretized PDE More...
 
pure real function getcfl (this)
 
subroutine fargoadvectionx (this, Fluxes, Mesh, Physics, Sources)
 Calculates the linear transport step in Fargo Scheme along x-axis [34] . More...
 
subroutine fargoadvectiony (this, Fluxes, Mesh, Physics, Sources)
 Calculates the linear transport step in Fargo Scheme along y-axis [34] . More...
 
subroutine fargoadvectionz (this, Fluxes, Mesh, Physics, Sources)
 Calculates the linear transport step in Fargo Scheme along z-axis [34] . More...
 
subroutine calcbackgroundvelocity (this, Mesh, Physics, pvar, cvar)
 Calculates new background velocity for fargo advection. More...
 
real function, dimension(mesh%igmin:mesh%igmax, mesh%jgmin:mesh%jgmax, mesh%kgmin:mesh%kgmax, physics%vdim) getcentrifugalvelocity (this, Mesh, Physics, Fluxes, Sources, dir_omega_, accel_, centrot)
 Compute velocity leading to a centrifugal acceleration with respect to some given axis of rotation which balances any local radial acceleration. If the acceleration is not given explicitly determine it from a complete evaluation of of the right hand side of the transport problem. More...
 
subroutine finalize_base (this)
 

Detailed Description

Author
Tobias Illenseer
Björn Sperling
Manuel Jung
Jannes Klee

Function/Subroutine Documentation

◆ acceptsolution()

subroutine timedisc_base_mod::acceptsolution ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources,
class(fluxes_base), intent(inout)  Fluxes,
real, intent(inout)  time,
real, intent(inout)  dt,
integer  iter 
)

Definition at line 911 of file timedisc_base.f90.

◆ adjusttimestep()

real function timedisc_base_mod::adjusttimestep ( class(timedisc_base), intent(inout)  this,
real, intent(in)  maxerr,
real, intent(in)  dtold 
)

adjust the time step

This function implements adaptive step size control based on an error estimate.

see E. Hairer, Solving Ordinary Differential Equ. II, 2ed, Springer (2.43c)

Definition at line 759 of file timedisc_base.f90.

Here is the call graph for this function:

◆ calcbackgroundvelocity()

subroutine timedisc_base_mod::calcbackgroundvelocity ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(marray_compound), intent(inout)  pvar,
class(marray_compound), intent(inout)  cvar 
)

Calculates new background velocity for fargo advection.

Attention
Only works when velocity is shifted in second direction.

Definition at line 2243 of file timedisc_base.f90.

◆ calctimestep()

real function timedisc_base_mod::calctimestep ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources,
class(fluxes_base), intent(inout)  Fluxes,
real, intent(in)  time,
integer, intent(inout)  dtcause 
)

Determines the CFL time step and time step limits due to source terms.

Definition at line 818 of file timedisc_base.f90.

◆ checkdata()

subroutine timedisc_base_mod::checkdata ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(fluxes_base), intent(in)  Fluxes,
class(marray_compound), intent(inout)  pvar,
class(marray_compound), intent(inout)  cvar,
integer, intent(in)  checkdatabm 
)

compute the RHS of the spatially discretized PDE

Definition at line 1666 of file timedisc_base.f90.

◆ computeerror()

real function timedisc_base_mod::computeerror ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(in)  Physics,
class(marray_compound), intent(inout)  cvar_high,
class(marray_compound), intent(inout)  cvar_low 
)
private

Definition at line 1012 of file timedisc_base.f90.

◆ computerhs()

subroutine timedisc_base_mod::computerhs ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources,
class(fluxes_base), intent(inout)  Fluxes,
real, intent(in)  time,
real, intent(in)  dt,
class(marray_compound), intent(inout)  pvar,
class(marray_compound), intent(inout)  cvar,
integer, intent(in)  checkdatabm,
class(marray_compound), intent(inout)  rhs 
)

compute the RHS of the spatially discretized PDE

This subroutine updates all data on the right hand side (RHS) of the ordinary differential equation obtained via spatial discretization of the system of partial differential equations (PDE). It expects that the conservative variables (cvar) contain valid data inside the computational domain. On the basis of this data the following update steps are performed in the given order:

  1. update of ghost cell data (this implies an update of the primitive variables (pvar) on the whole grid including ghost cells)
  2. update of intercell numerical fluxes (this implies an update of all face states and additional data necessary to compute the fluxes)
  3. update geometrical source terms
  4. update external source terms (this implies an update of all auxiliary data arrays used for sources terms)
Todo:
Very hacky implementation of pluto style angular momentum conservation. A better implementation is needed, but we would need to restucture the timedisc module

Definition at line 1073 of file timedisc_base.f90.

◆ fargoadvectionx()

subroutine timedisc_base_mod::fargoadvectionx ( class(timedisc_base), intent(inout)  this,
class(fluxes_base), intent(inout)  Fluxes,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources 
)

Calculates the linear transport step in Fargo Scheme along x-axis [34] .

The linear transport calculation is done by operator splitting. This yields an equation

\[ \frac{\partial q}{\partial t} + \mathbf{w}\cdot \nabla q = 0, \]

which can be descritized as

\[ q_j^{n+1}=q_j^{n}+ \frac{1}{\Delta y} \left( \int_{y_{j+\frac{1}{2}}-w\Delta t}^{y_{j+\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi - \int_{y_{j-\frac{1}{2}}-w\Delta t}^{y_{j-\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi \right). \]

FargoAdvectionX does the advection step along the x-axis with \( \mathbf{w} = w \mathbf{e_x} \).

Definition at line 1781 of file timedisc_base.f90.

◆ fargoadvectiony()

subroutine timedisc_base_mod::fargoadvectiony ( class(timedisc_base), intent(inout)  this,
class(fluxes_base), intent(inout)  Fluxes,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources 
)

Calculates the linear transport step in Fargo Scheme along y-axis [34] .

The linear transport calculation is done by operator splitting. This yields an equation

\[ \frac{\partial q}{\partial t} + \mathbf{w}\cdot \nabla q = 0, \]

which can be descritized as

\[ q_j^{n+1}=q_j^{n}+ \frac{1}{\Delta y} \left( \int_{y_{j+\frac{1}{2}}-w\Delta t}^{y_{j+\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi - \int_{y_{j-\frac{1}{2}}-w\Delta t}^{y_{j-\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi \right). \]

FargoAdvectionY does the advection step along the x-axis with \( \mathbf{w} = w \mathbf{e_y} \).

Definition at line 1935 of file timedisc_base.f90.

◆ fargoadvectionz()

subroutine timedisc_base_mod::fargoadvectionz ( class(timedisc_base), intent(inout)  this,
class(fluxes_base), intent(inout)  Fluxes,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources 
)

Calculates the linear transport step in Fargo Scheme along z-axis [34] .

The linear transport calculation is done by operator splitting. This yields an equation

\[ \frac{\partial q}{\partial t} + \mathbf{w}\cdot \nabla q = 0, \]

which can be descritized as

\[ q_j^{n+1}=q_j^{n}+ \frac{1}{\Delta y} \left( \int_{y_{j+\frac{1}{2}}-w\Delta t}^{y_{j+\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi - \int_{y_{j-\frac{1}{2}}-w\Delta t}^{y_{j-\frac{1}{2}}} q_j^{n}(\xi)\mathrm{d}\xi \right). \]

FargoAdvectionY does the advection step along the x-axis with \( \mathbf{w} = w \mathbf{e_y} \).

Definition at line 2089 of file timedisc_base.f90.

◆ finalize_base()

subroutine timedisc_base_mod::finalize_base ( class(timedisc_base this)
private

Definition at line 2532 of file timedisc_base.f90.

◆ getcentrifugalvelocity()

real function, dimension(mesh%igmin:mesh%igmax,mesh%jgmin:mesh%jgmax,mesh%kgmin:mesh%kgmax,physics%vdim) timedisc_base_mod::getcentrifugalvelocity ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(fluxes_base), intent(inout)  Fluxes,
class(sources_list), intent(inout), allocatable  Sources,
real, dimension(3), intent(in), optional  dir_omega_,
real, dimension(mesh%igmin:mesh%igmax,mesh%jgmin:mesh%jgmax,mesh%kgmin:mesh%kgmax,physics%vdim), intent(in), optional  accel_,
real, dimension(3), intent(in), optional  centrot 
)
private

Compute velocity leading to a centrifugal acceleration with respect to some given axis of rotation which balances any local radial acceleration. If the acceleration is not given explicitly determine it from a complete evaluation of of the right hand side of the transport problem.

The algorithm follows these steps:

  1. Get curvilinear components of the vectors pointing from the origin/center of rotation into each cell.
  2. Get curvilinear components of the vectors pointing from the axis of rotation perpendicular to that axis into each cell, i. e.

    \[ \vec{R} = \vec{r} - (\hat{e}_\Omega\cdot\vec{r})\hat{e}_\Omega \]

    . and the components of the azimuthal unit vector in terms of the local orthnormal basis

    \[ \hat{e}_\varphi = \hat{e}_\Omega \times \vec{r} / r \]

  3. Determine the acceleration to balance.
  4. Compute absolute value of azimuthal velocity using the balance law

    \[ \Omega^2 R\hat{e}_R = -\vec{a} \Leftrightarrow v_\varphi = \sqrt{\Omega^2 R^2} = \sqrt{-R \hat{e}_R\cdot\vec{a}} \]

  5. Compute components of \( v_\varphi \hat{e}_\varphi \) with respect to the local orthonormal basis.
Todo:
move code depending on physics into appropriate subroutines in physics modules

Definition at line 2355 of file timedisc_base.f90.

Here is the call graph for this function:

◆ getcfl()

pure real function timedisc_base_mod::getcfl ( class(timedisc_base), intent(in)  this)
private

Definition at line 1753 of file timedisc_base.f90.

◆ getorder()

pure integer function timedisc_base_mod::getorder ( class(timedisc_base), intent(in)  this)

Definition at line 1744 of file timedisc_base.f90.

Here is the caller graph for this function:

◆ inittimedisc()

subroutine timedisc_base_mod::inittimedisc ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(inout)  Mesh,
class(physics_base), intent(in)  Physics,
type(dict_typ), pointer  config,
type(dict_typ), pointer  IO,
integer, intent(in)  ttype,
character(len=32), intent(in)  tname 
)
private

Definition at line 224 of file timedisc_base.f90.

◆ integrationstep()

subroutine timedisc_base_mod::integrationstep ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(inout)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources,
class(fluxes_base), intent(inout)  Fluxes,
integer, intent(inout)  iter,
type(dict_typ), pointer  config,
type(dict_typ), pointer  IO 
)
private

Definition at line 649 of file timedisc_base.f90.

◆ rejectsolution()

subroutine timedisc_base_mod::rejectsolution ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(inout)  Physics,
class(sources_list), intent(inout), allocatable  Sources,
class(fluxes_base), intent(inout)  Fluxes,
real, intent(in)  time,
real, intent(inout)  dt 
)
private

Definition at line 961 of file timedisc_base.f90.

◆ setoutput()

subroutine timedisc_base_mod::setoutput ( class(timedisc_base), intent(inout)  this,
class(mesh_base), intent(in)  Mesh,
class(physics_base), intent(in)  Physics,
type(dict_typ), pointer  config,
type(dict_typ), pointer  IO 
)

Definition at line 534 of file timedisc_base.f90.

Variable Documentation

◆ cash_karp

integer, parameter, public timedisc_base_mod::cash_karp = 3

Definition at line 189 of file timedisc_base.f90.

◆ check_all

integer, parameter, public timedisc_base_mod::check_all = NOT(CHECK_NOTHING)

Definition at line 199 of file timedisc_base.f90.

◆ check_csound

integer, parameter, public timedisc_base_mod::check_csound = INT(b'1')

Definition at line 200 of file timedisc_base.f90.

◆ check_invalid

integer, parameter, public timedisc_base_mod::check_invalid = INT(b'1000')

Definition at line 203 of file timedisc_base.f90.

◆ check_nothing

integer, parameter, public timedisc_base_mod::check_nothing = INT(b'0')

Definition at line 198 of file timedisc_base.f90.

◆ check_pmin

integer, parameter, public timedisc_base_mod::check_pmin = INT(b'10')

Definition at line 201 of file timedisc_base.f90.

◆ check_rhomin

integer, parameter, public timedisc_base_mod::check_rhomin = INT(b'100')

Definition at line 202 of file timedisc_base.f90.

◆ check_tmin

integer, parameter, public timedisc_base_mod::check_tmin = INT(b'10000')

Definition at line 204 of file timedisc_base.f90.

◆ dormand_prince

integer, parameter, public timedisc_base_mod::dormand_prince = 4

Definition at line 190 of file timedisc_base.f90.

◆ dtcause_cfl

integer, parameter, public timedisc_base_mod::dtcause_cfl = 0

Definition at line 193 of file timedisc_base.f90.

◆ dtcause_erradj

integer, parameter, public timedisc_base_mod::dtcause_erradj = -1

Definition at line 194 of file timedisc_base.f90.

◆ dtcause_fileio

integer, parameter, public timedisc_base_mod::dtcause_fileio = -4

smallest ts due to fileio

Definition at line 196 of file timedisc_base.f90.

◆ dtcause_smallerr

integer, parameter, public timedisc_base_mod::dtcause_smallerr = -2

Definition at line 195 of file timedisc_base.f90.

◆ fargo_method

character(len=40), dimension(3), parameter timedisc_base_mod::fargo_method = (/ "dynamic background velocity ", "user supplied fixed background velocity ", "shearing box " /)
private

Definition at line 206 of file timedisc_base.f90.

◆ modified_euler

integer, parameter, public timedisc_base_mod::modified_euler = 1

Definition at line 187 of file timedisc_base.f90.

◆ rk_fehlberg

integer, parameter, public timedisc_base_mod::rk_fehlberg = 2

Definition at line 188 of file timedisc_base.f90.

◆ ssprk

integer, parameter, public timedisc_base_mod::ssprk = 5

Definition at line 191 of file timedisc_base.f90.