Monotone Lax-Friedrichs flux. It is more dissipative than the Godunov method, but computationally less demanding. Source: Equation 2.72, page 21.
Note
Although it might be useful, this procedure cannot be defined as elemental, because it has a dummy procedure as argument.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(flux) | :: | f |
flux function, |
|||
real(kind=rk), | intent(in) | :: | vm |
left (minus) reconstruction, |
||
real(kind=rk), | intent(in) | :: | vp |
right (plus) reconstruction, |
||
real(kind=rk), | intent(in) | :: | x(:) |
at flux interface, |
||
real(kind=rk), | intent(in) | :: | t |
time, |
||
real(kind=rk), | intent(in) | :: | alpha |
in the domain on the problem |
pure real(rk) function lax_friedrichs(f, vm, vp, x, t, alpha) result(res) !! Monotone Lax-Friedrichs flux. It is more dissipative than the Godunov method, but !! computationally less demanding. !! Source: Equation 2.72, page 21. !! !! @note !! Although it might be useful, this procedure cannot be defined as *elemental*, !! because it has a dummy procedure as argument. procedure(flux) :: f !! flux function, \( f(v, x, t) \) real(rk), intent(in) :: vm !! left (minus) reconstruction, \( v_{i+1/2}^- \) real(rk), intent(in) :: vp !! right (plus) reconstruction, \( v_{i+1/2}^+ = v_{(i+1)-1/2}^- \) real(rk), intent(in) :: x(:) !! \(x\) at flux interface, \( x_{i+1/2} \) real(rk), intent(in) :: t !! time, \( t \) real(rk), intent(in) :: alpha !! \( \max(|f'(v)|) \) in the domain on the problem res = (f(vm, x, t) + f(vp, x, t) - alpha*(vp - vm))/2 end function lax_friedrichs