EVOLUTION, FORCING: Add injected flag to update_forcing().

This flag indicates whether the energy injection was performed. If not,
there is no need for the variable update.

Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
Grzegorz Kowal 2023-08-29 17:43:03 -03:00
parent 55b13fb895
commit 26022c4ce0
2 changed files with 45 additions and 20 deletions

View File

@ -1097,6 +1097,8 @@ module evolution
integer, intent(out) :: status
logical :: injected
!-------------------------------------------------------------------------------
!
! advance the solution using the selected method
@ -1106,11 +1108,13 @@ module evolution
! add forcing contribution, requires the boundary and primitive variable update
!
if (forcing_enabled) then
call update_forcing(time, dt)
call update_forcing(time, dt, injected)
if (injected) then
call update_variables(time + dt, 0.0d+00, status)
if (status /= 0) go to 100
end if
end if
! check if we need to perform the refinement step
!

View File

@ -35,8 +35,9 @@ module forcing
! interfaces for procedure pointers
!
abstract interface
subroutine update_forcing_iface(t, dt)
real(kind=8), intent(in) :: t, dt
subroutine update_forcing_iface(t, dt, injected)
real(kind=8), intent(in ) :: t, dt
logical , intent( out) :: injected
end subroutine
end interface
@ -969,23 +970,25 @@ module forcing
! subroutine UPDATE_FORCING_EDDIES:
! --------------------------------
!
! Subroutine adds the energy injection terms.
! This subroutine adds the energy injection terms using simple eddy mixing.
!
! Arguments:
!
! t - time;
! dt - time step;
! injected - flag indicating whether the injection was performed;
!
!===============================================================================
!
subroutine update_forcing_eddies(t, dt)
subroutine update_forcing_eddies(t, dt, injected)
use coordinates, only : xmin, ymin, zmin, xlen, ylen, zlen
use random , only : randuni, randsym
implicit none
real(kind=8), intent(in) :: t, dt
real(kind=8), intent(in ) :: t, dt
logical , intent( out) :: injected
integer :: ni, n
#if NDIMS == 3
@ -995,6 +998,8 @@ module forcing
!-------------------------------------------------------------------------------
!
injected = .false.
if (t < tinj_start .or. t > tinj_stop) then
arms = 0.0d+00
dinj = 0.0d+00
@ -1042,6 +1047,8 @@ module forcing
end do
injected = .true.
end if
!-------------------------------------------------------------------------------
@ -1053,16 +1060,18 @@ module forcing
! subroutine UPDATE_FORCING_OH:
! ----------------------------
!
! Subroutine adds the energy injection terms for OrnsteinUhlenbeck driving.
! This subroutine adds the energy injection terms using OrnsteinUhlenbeck
! driving.
!
! Arguments:
!
! t - time;
! dt - time step;
! injected - flag indicating whether the injection was performed;
!
!===============================================================================
!
subroutine update_forcing_oh(t, dt)
subroutine update_forcing_oh(t, dt, injected)
#if NDIMS == 3
use constants, only : pi2
@ -1071,7 +1080,8 @@ module forcing
implicit none
real(kind=8), intent(in) :: t, dt
real(kind=8), intent(in ) :: t, dt
logical , intent( out) :: injected
integer :: l
real(kind=8) :: acoeff, dcoeff
@ -1084,6 +1094,8 @@ module forcing
!-------------------------------------------------------------------------------
!
injected = .false.
if (t < tinj_start .or. t > tinj_stop) then
arms = 0.0d+00
dinj = 0.0d+00
@ -1136,6 +1148,8 @@ module forcing
!
rinj = dinj / dt
injected = .true.
!-------------------------------------------------------------------------------
!
end subroutine update_forcing_oh
@ -1145,23 +1159,26 @@ module forcing
! subroutine UPDATE_FORCING_ALVELIUS:
! ----------------------------------
!
! Subroutine adds the energy injection terms for Alvelius driving.
! This subroutine adds the energy injection terms using Alvelius
! driving.
!
! Arguments:
!
! t - time;
! dt - time step;
! injected - flag indicating whether the injection was performed;
!
!===============================================================================
!
subroutine update_forcing_alvelius(t, dt)
subroutine update_forcing_alvelius(t, dt, injected)
use constants, only : pi2
use random , only : randuni
implicit none
real(kind=8), intent(in) :: t, dt
real(kind=8), intent(in ) :: t, dt
logical , intent( out) :: injected
integer :: l
real(kind=8) :: th1, dinj, sqdt
@ -1176,6 +1193,8 @@ module forcing
!-------------------------------------------------------------------------------
!
injected = .false.
if (t < tinj_start .or. t > tinj_stop) then
arms = 0.0d+00
dinj = 0.0d+00
@ -1258,6 +1277,8 @@ module forcing
!
rinj = dinj / dt
injected = .true.
!-------------------------------------------------------------------------------
!
end subroutine update_forcing_alvelius