timestep_module Module Reference

The timestep_module contains all timestep subroutines and newton-raphsons to solve the ODE system. More...

Functions/Subroutines

subroutine timestep (ctime, temp, dens, entropy, Y, dYdt, Ye, delYe, h, evmode)
 Network timestepping along the thermodynamic trajectory, specified by {ztime, ztemp, zdens}[zsteps]. More...
 
subroutine, private restrict_timestep (ctime, h, temp, dens, ttemp, tdens)
 Restricting the previously selected timestep. More...
 
real(r_kind) function, private ye_timestep (h, Ye, delYe)
 Estimate the timestep from the change of electron fraction. More...
 
real(r_kind) function, private abchange_timestep (h, Y, dYdt)
 Estimate the timestep from the abundances change. More...
 
real(r_kind) function, private hydro_timestep (h1, ctime, temp, dens, entropy, ttemp, tdens, Ye)
 Estimate the hydro timestep. More...
 
subroutine, public select_optimistic_timestep (init)
 Selects the timestep. More...
 
subroutine, private update_hydro (time_i, temperature, density, radius, entropy, efraction, T_tr)
 Returns temperature, density, and radius for a given time "time_i". More...
 
subroutine, public advance_implicit_euler (cnt)
 Advance system to the next step for the implicit Euler method. More...
 
subroutine, public advance_gear (cnt)
 Advance system to the next step using the Gear method (using gear_module) More...
 

Detailed Description

The timestep_module contains all timestep subroutines and newton-raphsons to solve the ODE system.

This module contains several different treatments for the timestep and for solving the ODE with different solvers (so far implicit euler and gear). The gear solver develops its own timestep in gear_module. It is only affected by the subroutine restrict_timestep.

Function/Subroutine Documentation

◆ abchange_timestep()

real(r_kind) function, private timestep_module::abchange_timestep ( real(r_kind), intent(in)  h,
real(r_kind), dimension(net_size), intent(in)  Y,
real(r_kind), dimension(net_size), intent(in)  dYdt 
)
private

Estimate the timestep from the abundances change.

In the network evolution mode, the timestep is calculated based on the most changing isotope. The timestep is the maximum of :

\[ h = \eta \cdot \left|\frac{Y}{\Delta Y} \right| \]

for all nuclei above a certain threshold abundance parameter_class::timestep_Ymin. The factor \( \eta \) is determined by the parameter parameter_class::timestep_factor.

Created: OK 28.08.2017

Parameters
[in]yabundances
[in]dydtchange of abundances
[in]htimescale on which this change happens

Definition at line 364 of file timestep_module.f90.

◆ advance_gear()

subroutine, public timestep_module::advance_gear ( integer, intent(in)  cnt)

Advance system to the next step using the Gear method (using gear_module)

This subroutine is similar to advance_implicit_euler. However, no additional check for convergence is done after the newton-raphson. This is not necessary as the Gear solver has its estimate of a possible error in the next step, given by evolving higher orders and comparing to lower ones.

See also
gear_module, Martin, D. 2017

Edited:

  • DM 15.06.2017
  • OK 27.08.2017
  • MR 22.12.2020
Parameters
[in]cntglobal iteration counter

Definition at line 819 of file timestep_module.f90.

Here is the call graph for this function:

◆ advance_implicit_euler()

subroutine, public timestep_module::advance_implicit_euler ( integer, intent(in)  cnt)

Advance system to the next step for the implicit Euler method.

This subroutine solves one timestep for the implicit Euler method via a newton-raphson (NR) scheme. The maximum amount of iterations for the NR is controlled by the parameter parameter_class::nr_maxcount. The NR is successful if the mass conservation deviates less than parameter_class::nr_tol. If the maximum iterations have been reached without convergence, it repeats the calculation with half of the timestep. This is done parameter_class::adapt_stepsize_maxcount many times. If after this amount of iterations still no convergence could be achieved it will raise an error.

See also
Winteler 2013

Edited:

  • OK: 15.06.2017 - Made into a separate subroutine
  • MR: 06.07.2022 - Let it crash if the result does not converge properly
Parameters
[in]cntglobal iteration counter

Definition at line 627 of file timestep_module.f90.

Here is the call graph for this function:

◆ hydro_timestep()

real(r_kind) function, private timestep_module::hydro_timestep ( real(r_kind), intent(in)  h1,
real(r_kind), intent(in)  ctime,
real(r_kind), intent(in)  temp,
real(r_kind), intent(in)  dens,
real(r_kind), intent(in)  entropy,
real(r_kind), intent(inout)  ttemp,
real(r_kind), intent(inout)  tdens,
real(r_kind), intent(in)  Ye 
)
private

Estimate the hydro timestep.

The timestep is restricted by the change of the hydrodynamic conditions (i.e., temperature and density). The maximum change of them is determined by the parameter parameter_class::timestep_hydro_factor. If the temperature or density change to rapidly the timestep will be decreased and ultimately, if the timestep is lower than \( 10^{-24}s \) an exception will be raised and the code terminates. The subroutine is also responsible to update the next temperature (ttemp) and density (tdens) that are used in restrict_timestep.

Note
This subroutine restricts only the timestep for the implicit euler solver, not for the gear solver.

Created: OK 28.08.2017

Parameters
[in]h1timescale on which this change happens
[in]ctimecurrent time
[in]temptemperature at ctime
[in]densdensity at ctime
[in]entropy[kB/baryon]
[in,out]ttempcomputed temperature at ctime + h
[in,out]tdenscomputed density at ctime + h
[in]yeelectron fraction

Definition at line 406 of file timestep_module.f90.

Here is the call graph for this function:

◆ restrict_timestep()

subroutine, private timestep_module::restrict_timestep ( real(r_kind), intent(in)  ctime,
real(r_kind), intent(inout)  h,
real(r_kind), intent(in)  temp,
real(r_kind), intent(in)  dens,
real(r_kind), intent(in)  ttemp,
real(r_kind), intent(in)  tdens 
)
private

Restricting the previously selected timestep.

Possible restrictions are given by the termination criterion and the trajectory datapoints as well as custom selected snapshots. Furthermore, it is restricted to hit the NSE transition temperature (parameter_class::nsetemp_cold).

Note
The timestep of the gear solver, which is evolved independently in gear_module is also effected by this subroutine.

Edited:

  • MR: 27.06.18
  • MR: 22.12.20 - Removed "const" parameter
Parameters
[in]ctimecurrent time
[in]temptemperature [GK]
[in]densdensity [g/cm^3]
[in]ttemptemperature [GK] of the next step
[in]tdensdensity [g/cm^3] of the next step
[in,out]hresulting timestep

Definition at line 140 of file timestep_module.f90.

Here is the call graph for this function:

◆ select_optimistic_timestep()

subroutine, public timestep_module::select_optimistic_timestep ( logical, intent(inout)  init)

Selects the timestep.

This subroutine is the main interface to the driver.f90 . It calls timestep and initializes the timestep in case of the first call.

Edited:

  • 17.12.2016

Definition at line 488 of file timestep_module.f90.

Here is the call graph for this function:

◆ timestep()

subroutine timestep_module::timestep ( real(r_kind), intent(in)  ctime,
real(r_kind), intent(in)  temp,
real(r_kind), intent(in)  dens,
real(r_kind), intent(in)  entropy,
real(r_kind), dimension(net_size), intent(in)  Y,
real(r_kind), dimension(net_size), intent(in)  dYdt,
real(r_kind), intent(in)  Ye,
real(r_kind), intent(in)  delYe,
real(r_kind), intent(inout)  h,
integer, intent(in)  evmode 
)

Network timestepping along the thermodynamic trajectory, specified by {ztime, ztemp, zdens}[zsteps].

This subroutine calls more specific subroutines to calculate the timestep. This also depends on the conditions. For example, in NSE it call ye_timestep, while it calls abchange_timestep otherwise.

Note
Gear solver will only be affected by the last call of this subroutine of restrict_timestep

Edited:

  • OK: 13.01.14
  • OK: 7.11.16
  • OK: 27.08.17 - merged in the rest of timestepping subroutines
  • MR: 27.06.18 - introduced additional subroutine "restrict_timestep"
Parameters
[in]ctimecurrent time
[in]temptemperature [GK]
[in]densdensity [g/cm^3]
[in]entropy[kB/baryon]
[in]yabundances
[in]dydtchange of abundances
[in]yeelectron fraction
[in]delyechange of elec.fraction (Ye-Ye_p)
[in,out]hresulting timestep
[in]evmodezone evolution mode

Definition at line 53 of file timestep_module.f90.

Here is the call graph for this function:

◆ update_hydro()

subroutine, private timestep_module::update_hydro ( real(r_kind), intent(in)  time_i,
real(r_kind), intent(inout)  temperature,
real(r_kind), intent(inout)  density,
real(r_kind), intent(inout)  radius,
real(r_kind), intent(in)  entropy,
real(r_kind), intent(in)  efraction,
real(r_kind), intent(out), optional  T_tr 
)
private

Returns temperature, density, and radius for a given time "time_i".

The time, entropy, and electron fraction are only inputs and are not changed. The update depends on the trajectory mode of the run and can be either interpolated from a trajectory, got from the expansion, or from an analytic expression.

Edited:

  • MR: 15.01.2020 - Moved it from advance_implicit_euler/gear to this new subroutine
  • MR: 10.02.2023 - Gave an optional output for the temperature
Parameters
[in]time_iTime [s]
[in]entropyEntropy [kB/nuc]
[in]efractionElectron fraction
[in,out]temperaturetemperature at time_i
[in,out]densitydensity at time_i
[in,out]radiusradius at time_i
[out]t_trtemperature at time_i, ignoring heating

Definition at line 530 of file timestep_module.f90.

Here is the call graph for this function:

◆ ye_timestep()

real(r_kind) function, private timestep_module::ye_timestep ( real(r_kind), intent(in)  h,
real(r_kind), intent(in)  Ye,
real(r_kind), intent(in)  delYe 
)
private

Estimate the timestep from the change of electron fraction.

In NSE, only (slow) weak reactions are taken into account when calculating the timestep. Therefore, only the change in \( y_e \) is considered here and the timestep is usually larger than in the full network mode. The timestep is calculated by:

\[ h = h_\mathrm{old} \cdot \eta \frac{y_e}{|\Delta y_e|} \]

with \( \eta \) being the parameter parameter_class::timestep_factor.

Created: OK 28.08.2017

Parameters
[in]yeelectron fraction
[in]delyechange of the electron fraction
[in]htimescale on which this change happens

Definition at line 335 of file timestep_module.f90.