Routines relative to time discretization of PDEs. More...
| Namespaces | |
| module | enum_flux_limiter | 
| Time order discrezation enumeration. | |
| module | enum_refined_velocity_method | 
| Time order discrezation enumeration. | |
| module | enum_time_order_discretization | 
| Time order discrezation enumeration. | |
| module | variables_time_discretization | 
| Declaration of variables associated to the time discretization. | |
| Functions | |
| subroutine | mod_compute_time_coefficients::compute_time_coefficients (time_order_discretization, time_step, time_step_n, time_discretization_coef_np1, time_discretization_coef_n, time_discretization_coef_nm1) | 
| Compute the time coefficients of the temporal scheme of a PDE. | |
| subroutine | compute_time_loop_cfl () | 
| Compute time steps that fulfill a multiple of the CFL. | |
| double precision function | mod_compute_time_step_cfl_advection::compute_time_step_cfl_advection (velocity, temporal_stability_factor) | 
| Compute the time step restricted by the CFL Warning: this is called "advection" refering the the mathematical advection equation Thus we look for the maximum velocity in the whole domain. | |
| double precision function | mod_compute_time_step_cfl_advection::compute_cfl_p_step_advection (velocity, time_step_local) | 
| Compute the time step restricted by the CFL for direction splitted advection Warning: this is called "advection" refering the the mathematical advection equation Thus we look for the maximum velocity in the whole domain. | |
| double precision function | mod_compute_time_step_cfl_advection::compute_cfl_acoustic (speed_of_sound, time_step_local) | 
| Compute the acoustic CFL. | |
| double precision function | mod_compute_time_step_cfl_advection::compute_time_step_cfl_acoustic (speed_of_sound, temporal_stability_factor) | 
| Compute the time_step related to acoustic CFL. | |
| double precision function | mod_compute_time_step_cfl_advection_surface_tension::compute_time_step_cfl_advection_surface_tension (velocity, density, temporal_stability_factor) | 
| Compute the time step restricted by the CFL Warning: this is called "advection" referring the the mathematical advection equation Thus we look for the maximum velocity plus the maximum capillary velocity (surface tension) in the whole domain. | |
| double precision function | mod_compute_time_step_cfl_advection_surface_tension::compute_time_step_cfl_advection_1d_surface_tension (velocity, density, temporal_stability_factor) | 
| Compute the time step restricted by the CFL Warning: this is called "advection" referring the the mathematical advection equation Thus we look for the maximum velocity in the whole domain. | |
| subroutine | stop_tests () | 
| Stops simulation when a stop tests is satisfied. | |
| subroutine | time_step_switch () | 
| Switch time dependant variables to next time step. | |
Routines relative to time discretization of PDEs.
The reader willl find in this directory routines to:
time_loop_cfl.f90)stop_tests.f90)time_step_switch.f90)The file variables.f90 contains the declaration of scalar variables associated to time discretization, and enum.f90 the enumarations. 
| double precision function mod_compute_time_step_cfl_advection::compute_cfl_acoustic | ( | double precision, dimension(:,:,:), intent(in), allocatable | speed_of_sound, | 
| double precision, intent(in) | time_step_local ) | 
Compute the acoustic CFL.
| [in] | speed_of_sound | the speed_of_sound field | 
| [in] | time_step_local | current time step | 
| double precision function mod_compute_time_step_cfl_advection::compute_cfl_p_step_advection | ( | type(t_face_field), intent(in) | velocity, | 
| double precision, intent(in) | time_step_local ) | 
Compute the time step restricted by the CFL for direction splitted advection Warning: this is called "advection" refering the the mathematical advection equation Thus we look for the maximum velocity in the whole domain.
| [in] | velocity | the velocity field | 
| [in] | time_step_local | current time step | 
| subroutine mod_compute_time_coefficients::compute_time_coefficients | ( | integer, intent(in) | time_order_discretization, | 
| double precision, intent(in) | time_step, | ||
| double precision, intent(in) | time_step_n, | ||
| double precision, intent(out) | time_discretization_coef_np1, | ||
| double precision, intent(out) | time_discretization_coef_n, | ||
| double precision, intent(out) | time_discretization_coef_nm1 ) | 
Compute the time coefficients of the temporal scheme of a PDE.
This temporal term is written as follow:
∂φ/∂t ≈ α φⁿ⁺¹ + β φⁿ + γ φⁿ⁻¹
1 -1α = ----— β = ----— γ = 0 Δtⁿ⁺¹ Δtⁿ⁺¹
2Δtⁿ⁺¹ + Δtⁿ -(Δtⁿ⁺¹ + Δtⁿ) Δtⁿ⁺¹α = ---------------— β = -------------— γ = -------------— Δtⁿ⁺¹(Δtⁿ⁺¹+Δtⁿ) Δtⁿ⁺¹Δtⁿ Δtⁿ(Δtⁿ⁺¹+Δtⁿ)
| subroutine compute_time_loop_cfl | 
Compute time steps that fulfill a multiple of the CFL.
This routine computes the time step that is equal to a multiple of the CFL.
This routine computes phase_advection_time_step, energy_time_step, navier_time_step equal to CFL*temporal_stability_factor 
| double precision function mod_compute_time_step_cfl_advection::compute_time_step_cfl_acoustic | ( | double precision, dimension(:,:,:), intent(in), allocatable | speed_of_sound, | 
| double precision, intent(in) | temporal_stability_factor ) | 
Compute the time_step related to acoustic CFL.
| [in] | speed_of_sound | the velocity field | 
| [in] | temporal_stability_factor | the coefficient multiplying the maximum time step | 
| double precision function mod_compute_time_step_cfl_advection::compute_time_step_cfl_advection | ( | type(t_face_field), intent(in) | velocity, | 
| double precision, intent(in) | temporal_stability_factor ) | 
Compute the time step restricted by the CFL Warning: this is called "advection" refering the the mathematical advection equation Thus we look for the maximum velocity in the whole domain.
| [in] | velocity | the velocity field | 
| [in] | temporal_stability_factor | the coefficient multiplying the maximum time step | 
| double precision function mod_compute_time_step_cfl_advection_surface_tension::compute_time_step_cfl_advection_1d_surface_tension | ( | type(t_face_field), intent(in) | velocity, | 
| double precision, dimension(:,:,:), intent(in) | density, | ||
| double precision, intent(in) | temporal_stability_factor ) | 
Compute the time step restricted by the CFL Warning: this is called "advection" referring the the mathematical advection equation Thus we look for the maximum velocity in the whole domain.
| [in] | velocity | the velocity field | 
| [in] | density | the density field | 
| [in] | temporal_stability_factor | the coefficient multiplying the maximum time step | 
| double precision function mod_compute_time_step_cfl_advection_surface_tension::compute_time_step_cfl_advection_surface_tension | ( | type(t_face_field), intent(in) | velocity, | 
| double precision, dimension(:,:,:), intent(in) | density, | ||
| double precision, intent(in) | temporal_stability_factor ) | 
Compute the time step restricted by the CFL Warning: this is called "advection" referring the the mathematical advection equation Thus we look for the maximum velocity plus the maximum capillary velocity (surface tension) in the whole domain.
| [in] | velocity | the velocity field | 
| [in] | density | the density field | 
| [in] | temporal_stability_factor | the coefficient multiplying the maximum time step | 
| subroutine stop_tests | 
Stops simulation when a stop tests is satisfied.
The stop tests are based on:
When all stop tests are satisfied, is_stop_simulation is set to .true. and the simulation is ended. 
| subroutine time_step_switch | 
Switch time dependant variables to next time step.
This routine switches time dependant variable defined at \(t^{n}\) to their value at \( t^{n+1} \). If second order time is used, switch is done between \(t^{n}\) and \(t^{n-1}\).