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

View File

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