ROMS
Loading...
Searching...
No Matches
bvf_mix.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#if defined NONLINEAR && defined BVF_MIXING
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This routine computes tracer vertical mixing as a function of the !
13! Brunt-Vaisala frequency. The vertical mixing of momentum is set !
14! to its background value. If static unstable regime, the vertical !
15! mixing is set to "bvf_nu0c". !
16! !
17!======================================================================!
18!
19 implicit none
20!
21 PRIVATE
22 PUBLIC :: bvf_mix
23!
24 CONTAINS
25!
26!***********************************************************************
27 SUBROUTINE bvf_mix (ng, tile)
28!***********************************************************************
29!
30 USE mod_param
31 USE mod_mixing
32!
33! Imported variable declarations.
34!
35 integer, intent(in) :: ng, tile
36!
37! Local variable declarations.
38!
39# include "tile.h"
40!
41 CALL bvf_mix_tile (ng, tile, &
42 & lbi, ubi, lbj, ubj, &
43 & imins, imaxs, jmins, jmaxs, &
44 & mixing(ng) % bvf, &
45 & mixing(ng) % AKt, &
46 & mixing(ng) % AKv)
47 RETURN
48 END SUBROUTINE bvf_mix
49!
50!***********************************************************************
51 SUBROUTINE bvf_mix_tile (ng, tile, &
52 & LBi, UBi, LBj, UBj, &
53 & IminS, ImaxS, JminS, JmaxS, &
54 & bvf, Akt, Akv)
55!***********************************************************************
56!
57 USE mod_param
58 USE mod_scalars
59!
61# ifdef DISTRIBUTE
63# endif
64!
65! Imported variable declarations.
66!
67 integer, intent(in) :: ng, tile
68 integer, intent(in) :: LBi, UBi, LBj, UBj
69 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
70!
71# ifdef ASSUMED_SHAPE
72 real(r8), intent(in) :: bvf(LBi:,LBj:,0:)
73 real(r8), intent(out) :: Akt(LBi:,LBj:,0:,:)
74 real(r8), intent(out) :: Akv(LBi:,LBj:,0:)
75# else
76 real(r8), intent(in) :: bvf(LBi:UBi,LBj:UBj,0:N(ng))
77 real(r8), intent(out) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
78 real(r8), intent(out) :: Akv(LBi:UBi,LBj:UBj,0:N(ng))
79# endif
80!
81! Local variable declarations.
82!
83 integer :: i, itrc, j, k
84
85 real(r8) :: cff
86
87# include "set_bounds.h"
88!
89!-----------------------------------------------------------------------
90! Set tracer diffusivity as function of the Brunt-vaisala frequency.
91! Set vertical viscosity to its background value.
92!-----------------------------------------------------------------------
93!
94 DO k=1,n(ng)-1
95 DO j=jstr,jend
96 DO i=istr,iend
97 akv(i,j,k)=akv_bak(ng)
98 IF (bvf(i,j,k).lt.0.0_r8) THEN
99 akv(i,j,k)=bvf_nu0c
100 akt(i,j,k,itemp)=bvf_nu0c
101# ifdef SALINITY
102 akt(i,j,k,isalt)=bvf_nu0c
103# endif
104 ELSE IF (bvf(i,j,k).eq.0.0_r8) THEN
105 akv(i,j,k)=akv_bak(ng)
106 akt(i,j,k,itemp)=akt_bak(itemp,ng)
107# ifdef SALINITY
108 akt(i,j,k,isalt)=akt_bak(isalt,ng)
109# endif
110 ELSE
111 cff=bvf_nu0/sqrt(bvf(i,j,k))
112 akt(i,j,k,itemp)=min(bvf_numax,max(bvf_numin,cff))
113 akv(i,j,k)=akt(i,j,k,itemp)
114# ifdef SALINITY
115 akt(i,j,k,isalt)=akt(i,j,k,itemp)
116# endif
117 END IF
118 END DO
119 END DO
120 END DO
121!
122!-----------------------------------------------------------------------
123! Exchange boundary data.
124!-----------------------------------------------------------------------
125!
126 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
127 CALL exchange_w3d_tile (ng, tile, &
128 & lbi, ubi, lbj, ubj, 0, n(ng), &
129 & akv)
130 DO itrc=1,nat
131 CALL exchange_w3d_tile (ng, tile, &
132 & lbi, ubi, lbj, ubj, 0, n(ng), &
133 & akt(:,:,:,itrc))
134 END DO
135 END IF
136
137# ifdef DISTRIBUTE
138 CALL mp_exchange3d (ng, tile, inlm, 1, &
139 & lbi, ubi, lbj, ubj, 0, n(ng), &
140 & nghostpoints, &
141 & ewperiodic(ng), nsperiodic(ng), &
142 & akv)
143 CALL mp_exchange4d (ng, tile, inlm, 1, &
144 & lbi, ubi, lbj, ubj, 0, n(ng), 1, nat, &
145 & nghostpoints, &
146 & ewperiodic(ng), nsperiodic(ng), &
147 & akt)
148# endif
149
150 RETURN
151 END SUBROUTINE bvf_mix_tile
152#endif
153 END MODULE bvf_mix_mod
subroutine bvf_mix_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, bvf, akt, akv)
Definition bvf_mix.F:55
subroutine, public bvf_mix(ng, tile)
Definition bvf_mix.F:28
subroutine exchange_w3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer, parameter inlm
Definition mod_param.F:662
integer nghostpoints
Definition mod_param.F:710
real(r8) bvf_numax
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(r8) bvf_numin
real(r8) bvf_nu0c
integer isalt
real(r8), dimension(:,:), allocatable akt_bak
real(r8), dimension(:), allocatable akv_bak
integer itemp
real(r8) bvf_nu0
subroutine mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, a, b, c)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)