diff --git a/README.md b/README.md index 44f2749..2e75a08 100644 --- a/README.md +++ b/README.md @@ -19,7 +19,7 @@ following features are already implemented: * high order reconstructions: from 2nd to 9th order WENO and MP, both explicit and compact methods, the 2nd order TVD interpolation has a number of limiters supported, -* Riemann solvers of Roe- and HLL-types (HLL, HLLC, and HLLD), +* Riemann solvers of KEPES-, Roe- and HLL-types (HLL, HLLC, and HLLD), * standard boundary conditions: periodic, open, reflective, hydrostatic, etc. * turbulence driving using Alvelius or Ornstein–Uhlenbeck methods, * viscous and resistive source terms, diff --git a/sources/forcing.F90 b/sources/forcing.F90 index e93e7c3..4916e49 100644 --- a/sources/forcing.F90 +++ b/sources/forcing.F90 @@ -96,7 +96,7 @@ module forcing real(kind=8), save :: fpower = 1.0d+00 real(kind=8), save :: lscale = 5.0d-01 real(kind=8), save :: tscale = 1.0d+00 - real(kind=8), save :: vscale = 5.0d-01 + real(kind=8), save :: vscale = 0.0d+00 real(kind=8), save :: ascale = 5.0d-01 real(kind=8), save :: kf = 2.0d+00 real(kind=8), save :: kl = 0.0d+00 @@ -274,6 +274,15 @@ module forcing end if tscale = 1.0d+00 end if + call get_parameter("driving_velocity", vscale) + if (vscale < 0.0d+00) then + if (verbose) then + write(*,*) + write(*,wfmt) "'driving_velocity' must be positive!" + write(*,wfmt) "Resetting it to the default value: 0.0!" + end if + vscale = 0.0d+00 + end if case("Alvelius", "alvelius", "al") forcing_enabled = .true. @@ -292,6 +301,15 @@ module forcing end if fpower = 1.0d+00 end if + call get_parameter("driving_velocity", vscale) + if (vscale < 0.0d+00) then + if (verbose) then + write(*,*) + write(*,wfmt) "'driving_velocity' must be positive!" + write(*,wfmt) "Resetting it to the default value: 0.0!" + end if + vscale = 0.0d+00 + end if case default injection_method = injection_none @@ -706,7 +724,11 @@ module forcing ! instantenous injection rate) ! lscale = 1.0d+00 / kf - vscale = lscale / tscale + if (vscale > 0.0d+00) then + tscale = lscale / vscale + else + vscale = lscale / tscale + end if ascale = sqrt(2.0d+00) * lscale / tscale**2 / sqrt(tscale) fpower = vscale**2 / tscale @@ -742,8 +764,13 @@ module forcing ! get the characteristic driving length, velocity, time, and acceleration scales ! lscale = 1.0d+00 / kf - tscale = sqrt(lscale / sqrt(2.0d+00 * fpower)) - vscale = lscale / tscale + if (vscale > 0.0d+00) then + tscale = lscale / vscale + fpower = 5.0d-01 * (vscale**2 / lscale)**2 + else + tscale = sqrt(lscale / sqrt(2.0d+00 * fpower)) + vscale = lscale / tscale + end if ascale = vscale / tscale ! normalize the driving amplitude profile to the driving power