ROMS
Loading...
Searching...
No Matches
lmd_vmix_mod Module Reference

Functions/Subroutines

subroutine, public lmd_vmix (ng, tile)
 
subroutine lmd_vmix_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nstp, hz, rho, u, v, t, alfaobeta, bvf, akt, akv)
 
subroutine lmd_finish (ng, tile)
 
subroutine lmd_finish_tile (ng, tile, lbi, ubi, lbj, ubj, nstp, hz, bvf, akt, akv)
 

Function/Subroutine Documentation

◆ lmd_finish()

subroutine lmd_vmix_mod::lmd_finish ( integer, intent(in) ng,
integer, intent(in) tile )
private

Definition at line 437 of file lmd_vmix.F.

438!***********************************************************************
439!
440 USE mod_param
441 USE mod_grid
442 USE mod_mixing
443 USE mod_ocean
444 USE mod_stepping
445!
446! Imported variable declarations.
447!
448 integer, intent(in) :: ng, tile
449!
450! Local variable declarations.
451!
452# include "tile.h"
453!
454 CALL lmd_finish_tile (ng, tile, &
455 & lbi, ubi, lbj, ubj, &
456 & nstp(ng), &
457 & grid(ng) % Hz, &
458 & mixing(ng) % bvf, &
459 & mixing(ng) % Akt, &
460 & mixing(ng) % Akv)
461 RETURN
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
integer, dimension(:), allocatable nstp

References mod_grid::grid, lmd_finish_tile(), mod_mixing::mixing, and mod_stepping::nstp.

Referenced by lmd_vmix().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ lmd_finish_tile()

subroutine lmd_vmix_mod::lmd_finish_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) nstp,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbi:,lbj:,0:), intent(in) bvf,
real(r8), dimension(lbi:,lbj:,0:,:), intent(inout) akt,
real(r8), dimension(lbi:,lbj:,0:), intent(inout) akv )
private

Definition at line 465 of file lmd_vmix.F.

470!***********************************************************************
471!
472 USE mod_param
473 USE mod_scalars
474!
475 USE bc_3d_mod, ONLY : bc_w3d_tile
476# ifdef DISTRIBUTE
478# endif
479!
480! Imported variable declarations.
481!
482 integer, intent(in) :: ng, tile
483 integer, intent(in) :: LBi, UBi, LBj, UBj
484 Integer, intent(in) :: nstp
485!
486# ifdef ASSUMED_SHAPE
487 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
488 real(r8), intent(in) :: bvf(LBi:,LBj:,0:)
489
490 real(r8), intent(inout) :: Akt(LBi:,LBj:,0:,:)
491 real(r8), intent(inout) :: Akv(LBi:,LBj:,0:)
492# else
493 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
494 real(r8), intent(in) :: bvf(LBi:UBi,LBj:UBj,0:N(ng))
495
496 real(r8), intent(inout) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
497 real(r8), intent(inout) :: Akv(LBi:UBi,LBj:UBj,0:N(ng))
498# endif
499!
500! Local variable declarations.
501!
502 integer :: i, itrc, j, k
503
504 real(r8), parameter :: eps = 1.0e-14_r8
505
506 real(r8) :: cff, lmd_iwm, lmd_iws, nu_sx, nu_sxc, shear2
507
508# include "set_bounds.h"
509!
510!-----------------------------------------------------------------------
511! Compute "interior" viscosities and diffusivities everywhere as
512! the superposition of three processes: local Richardson number
513! instability due to resolved vertical shear, internal wave
514! breaking, and double diffusion.
515!-----------------------------------------------------------------------
516!
517 DO k=1,n(ng)-1
518 DO j=jstr,jend
519 DO i=istr,iend
520!
521! Compute interior convective diffusivity due to static instability
522! mixing.
523!
524# ifdef LMD_CONVEC
525 cff=max(bvf(i,j,k),lmd_bvfcon)
526 cff=min(1.0_r8,(lmd_bvfcon-cff)/lmd_bvfcon)
527 nu_sxc=1.0_r8-cff*cff
528 nu_sxc=nu_sxc*nu_sxc*nu_sxc
529# else
530 nu_sxc=0.0_r8
531# endif
532!
533! Sum contributions due to internal wave breaking, shear instability
534! and convective diffusivity due to shear instability.
535!
536 akv(i,j,k)=akv(i,j,k)+lmd_nu0c*nu_sxc
537 akt(i,j,k,itemp)=akt(i,j,k,itemp)+lmd_nu0c*nu_sxc
538# ifdef SALINITY
539 akt(i,j,k,isalt)=akt(i,j,k,isalt)+lmd_nu0c*nu_sxc
540# endif
541# if defined LIMIT_VDIFF || defined LIMIT_VVISC
542!
543! Limit vertical mixing coefficients with the upper threshold value.
544! This is an engineering fix but it can be based on the fact that
545! vertical mixing in the ocean from indirect observations is not
546! higher than the threshold value.
547!
548# ifdef LIMIT_VVISC
549 akv(i,j,k)=min(akv_limit(ng), akv(i,j,k))
550# endif
551# ifdef LIMIT_VDIFF
552 akt(i,j,k,itemp)=min(akt_limit(itemp,ng), akt(i,j,k,itemp))
553# ifdef SALINITY
554 akt(i,j,k,isalt)=min(akt_limit(isalt,ng), akt(i,j,k,isalt))
555# endif
556# endif
557# endif
558 END DO
559 END DO
560 END DO
561!
562! Apply boundary conditions.
563!
564 DO k=0,n(ng)
565 IF (domain(ng)%Western_Edge(tile)) THEN
566 DO j=jstr,jend
567 DO itrc=1,nat
568 akt(istr-1,j,k,itrc)=akt(istr,j,k,itrc)
569 END DO
570 akv(istr-1,j,k)=akv(istr,j,k)
571 END DO
572 END IF
573 IF (domain(ng)%Eastern_Edge(tile)) THEN
574 DO j=jstr,jend
575 DO itrc=1,nat
576 akt(iend+1,j,k,itrc)=akt(iend,j,k,itrc)
577 END DO
578 akv(iend+1,j,k)=akv(iend,j,k)
579 END DO
580
581 END IF
582 IF (domain(ng)%Southern_Edge(tile)) THEN
583 DO i=istr,iend
584 DO itrc=1,nat
585 akt(i,jstr-1,k,itrc)=akt(i,jstr,k,itrc)
586 END DO
587 akv(i,jstr-1,k)=akv(i,jstr,k)
588 END DO
589 END IF
590 IF (domain(ng)%Northern_Edge(tile)) THEN
591 DO i=istr,iend
592 DO itrc=1,nat
593 akt(i,jend+1,k,itrc)=akt(i,jend,k,itrc)
594 END DO
595 akv(i,jend+1,k)=akv(i,jend,k)
596 END DO
597 END IF
598 IF (domain(ng)%SouthWest_Corner(tile)) THEN
599 DO itrc=1,nat
600 akt(istr-1,jstr-1,k,itrc)=0.5_r8* &
601 & (akt(istr ,jstr-1,k,itrc)+ &
602 & akt(istr-1,jstr ,k,itrc))
603 END DO
604 akv(istr-1,jstr-1,k)=0.5_r8* &
605 & (akv(istr ,jstr-1,k)+ &
606 & akv(istr-1,jstr ,k))
607 END IF
608 IF (domain(ng)%SouthEast_Corner(tile)) THEN
609 DO itrc=1,nat
610 akt(iend+1,jstr-1,k,itrc)=0.5_r8* &
611 & (akt(iend ,jstr-1,k,itrc)+ &
612 & akt(iend+1,jstr ,k,itrc))
613 END DO
614 akv(iend+1,jstr-1,k)=0.5_r8* &
615 & (akv(iend ,jstr-1,k)+ &
616 & akv(iend+1,jstr ,k))
617 END IF
618 IF (domain(ng)%NorthWest_Corner(tile)) THEN
619 DO itrc=1,nat
620 akt(istr-1,jend+1,k,itrc)=0.5_r8* &
621 & (akt(istr ,jend+1,k,itrc)+ &
622 & akt(istr-1,jend ,k,itrc))
623 END DO
624 akv(istr-1,jend+1,k)=0.5_r8* &
625 & (akv(istr ,jend+1,k)+ &
626 & akv(istr-1,jend ,k))
627 END IF
628 IF (domain(ng)%NorthEast_Corner(tile)) THEN
629 DO itrc=1,nat
630 akt(iend+1,jend+1,k,itrc)=0.5_r8* &
631 & (akt(iend ,jend+1,k,itrc)+ &
632 & akt(iend+1,jend ,k,itrc))
633 END DO
634 akv(iend+1,jend+1,k)=0.5_r8* &
635 & (akv(iend ,jend+1,k)+ &
636 & akv(iend+1,jend ,k))
637 END IF
638 END DO
639!
640 CALL bc_w3d_tile (ng, tile, &
641 & lbi, ubi, lbj, ubj, 0, n(ng), &
642 & akv)
643 DO itrc=1,nat
644 CALL bc_w3d_tile (ng, tile, &
645 & lbi, ubi, lbj, ubj, 0, n(ng), &
646 & akt(:,:,:,itrc))
647 END DO
648# ifdef DISTRIBUTE
649 CALL mp_exchange3d (ng, tile, inlm, 1, &
650 & lbi, ubi, lbj, ubj, 0, n(ng), &
651 & nghostpoints, &
652 & ewperiodic(ng), nsperiodic(ng), &
653 & akv)
654 CALL mp_exchange4d (ng, tile, inlm, 1, &
655 & lbi, ubi, lbj, ubj, 0, n(ng), 1, nat, &
656 & nghostpoints, &
657 & ewperiodic(ng), nsperiodic(ng), &
658 & akt)
659# endif
660!
661 RETURN
subroutine bc_w3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, a)
Definition bc_3d.F:591
integer nat
Definition mod_param.F:499
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
real(r8) lmd_nu0c
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(r8) lmd_bvfcon
integer isalt
integer itemp
real(r8), dimension(:,:), allocatable akt_limit
real(r8), dimension(:), allocatable akv_limit
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)

References mod_scalars::akt_limit, mod_scalars::akv_limit, bc_3d_mod::bc_w3d_tile(), mod_param::domain, mod_scalars::ewperiodic, mod_param::inlm, mod_scalars::isalt, mod_scalars::itemp, mod_scalars::lmd_bvfcon, mod_scalars::lmd_nu0c, mp_exchange_mod::mp_exchange3d(), mp_exchange_mod::mp_exchange4d(), mod_param::nghostpoints, and mod_scalars::nsperiodic.

Referenced by lmd_finish().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ lmd_vmix()

subroutine, public lmd_vmix_mod::lmd_vmix ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 33 of file lmd_vmix.F.

34!***********************************************************************
35!
36 USE mod_param
37 USE mod_grid
38 USE mod_mixing
39 USE mod_ocean
40 USE mod_scalars
41 USE mod_stepping
42!
43# ifdef LMD_SKPP
44 USE lmd_skpp_mod, ONLY : lmd_skpp
45# endif
46# ifdef LMD_BKPP
47 USE lmd_bkpp_mod, ONLY : lmd_bkpp
48# endif
49!
50! Imported variable declarations.
51!
52 integer, intent(in) :: ng, tile
53!
54! Local variable declarations.
55!
56 character (len=*), parameter :: MyFile = &
57 & __FILE__
58!
59# include "tile.h"
60!
61# ifdef PROFILE
62 CALL wclock_on (ng, inlm, 18, __line__, myfile)
63# endif
64 IF (.not.perfectrst(ng).or.iic(ng).ne.ntstart(ng)) THEN
65 CALL lmd_vmix_tile (ng, tile, &
66 & lbi, ubi, lbj, ubj, &
67 & imins, imaxs, jmins, jmaxs, &
68 & nstp(ng), &
69 & grid(ng) % Hz, &
70# ifndef RI_SPLINES
71 & grid(ng) % z_r, &
72# endif
73 & ocean(ng) % rho, &
74 & ocean(ng) % u, &
75 & ocean(ng) % v, &
76# ifdef LMD_DDMIX
77 & ocean(ng) % t, &
78 & mixing(ng) % alfaobeta, &
79# endif
80 & mixing(ng) % bvf, &
81 & mixing(ng) % Akt, &
82 & mixing(ng) % Akv)
83# ifdef LMD_SKPP
84 CALL lmd_skpp (ng, tile)
85# endif
86# ifdef LMD_BKPP
87 CALL lmd_bkpp (ng, tile)
88# endif
89 CALL lmd_finish (ng, tile)
90 END IF
91# ifdef PROFILE
92 CALL wclock_off (ng, inlm, 18, __line__, myfile)
93# endif
94!
95 RETURN
subroutine lmd_bkpp(ng, tile)
Definition lmd_bkpp.F:45
subroutine lmd_skpp(ng, tile)
Definition lmd_skpp.F:45
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable iic
logical, dimension(:), allocatable perfectrst
integer, dimension(:), allocatable ntstart
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_grid::grid, mod_scalars::iic, mod_param::inlm, lmd_bkpp(), lmd_finish(), lmd_skpp(), lmd_vmix_tile(), mod_mixing::mixing, mod_stepping::nstp, mod_scalars::ntstart, mod_ocean::ocean, mod_scalars::perfectrst, wclock_off(), and wclock_on().

Referenced by main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ lmd_vmix_tile()

subroutine lmd_vmix_mod::lmd_vmix_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nstp,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbi:,lbj:,:), intent(in) rho,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) u,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) v,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) t,
real(r8), dimension(lbi:,lbj:,0:), intent(in) alfaobeta,
real(r8), dimension(lbi:,lbj:,0:), intent(in) bvf,
real(r8), dimension(lbi:,lbj:,0:,:), intent(inout) akt,
real(r8), dimension(lbi:,lbj:,0:), intent(inout) akv )
private

Definition at line 99 of file lmd_vmix.F.

112!***********************************************************************
113!
114 USE mod_param
115 USE mod_scalars
116!
117! Imported variable declarations.
118!
119 integer, intent(in) :: ng, tile
120 integer, intent(in) :: LBi, UBi, LBj, UBj
121 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
122 Integer, intent(in) :: nstp
123!
124# ifdef ASSUMED_SHAPE
125 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
126# ifndef RI_SPLINES
127 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
128# endif
129 real(r8), intent(in) :: rho(LBi:,LBj:,:)
130 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
131 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
132# ifdef LMD_DDMIX
133 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
134 real(r8), intent(in) :: alfaobeta(LBi:,LBj:,0:)
135# endif
136 real(r8), intent(in) :: bvf(LBi:,LBj:,0:)
137
138 real(r8), intent(inout) :: Akt(LBi:,LBj:,0:,:)
139 real(r8), intent(inout) :: Akv(LBi:,LBj:,0:)
140# else
141 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
142# ifndef RI_SPLINES
143 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
144# endif
145 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
146 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),3)
147 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),3)
148# ifdef LMD_DDMIX
149 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
150 real(r8), intent(in) :: alfaobeta(LBi:UBi,LBj:UBj,0:N(ng))
151# endif
152 real(r8), intent(in) :: bvf(LBi:UBi,LBj:UBj,0:N(ng))
153
154 real(r8), intent(inout) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
155 real(r8), intent(inout) :: Akv(LBi:UBi,LBj:UBj,0:N(ng))
156# endif
157!
158! Local variable declarations.
159!
160 integer :: i, itrc, j, k
161
162 real(r8), parameter :: eps = 1.0e-14_r8
163
164 real(r8) :: cff, lmd_iwm, lmd_iws, nu_sx, nu_sxc, shear2
165# ifdef LMD_DDMIX
166 real(r8) :: Rrho, ddDS, ddDT, nu_dds, nu_ddt
167# endif
168
169 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: Rig
170
171 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC
172 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dR
173 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dU
174 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dV
175
176# include "set_bounds.h"
177
178# ifdef LMD_RIMIX
179!
180!-----------------------------------------------------------------------
181! Compute gradient Richardson number.
182!-----------------------------------------------------------------------
183!
184! Compute gradient Richardson number at horizontal RHO-points and
185! vertical W-points. If zero or very small velocity shear, bound
186! computation by a large negative value.
187!
188# ifdef RI_SPLINES
189 DO j=max(1,jstr-1),min(jend+1,mm(ng))
190 DO i=max(1,istr-1),min(iend+1,lm(ng))
191 fc(i,0)=0.0_r8
192 dr(i,0)=0.0_r8
193 du(i,0)=0.0_r8
194 dv(i,0)=0.0_r8
195 END DO
196 DO k=1,n(ng)-1
197 DO i=max(1,istr-1),min(iend+1,lm(ng))
198 cff=1.0_r8/(2.0_r8*hz(i,j,k+1)+ &
199 & hz(i,j,k)*(2.0_r8-fc(i,k-1)))
200 fc(i,k)=cff*hz(i,j,k+1)
201 dr(i,k)=cff*(6.0_r8*(rho(i,j,k+1)-rho(i,j,k))- &
202 & hz(i,j,k)*dr(i,k-1))
203 du(i,k)=cff*(3.0_r8*(u(i ,j,k+1,nstp)-u(i ,j,k,nstp)+ &
204 & u(i+1,j,k+1,nstp)-u(i+1,j,k,nstp))- &
205 & hz(i,j,k)*du(i,k-1))
206 dv(i,k)=cff*(3.0_r8*(v(i,j ,k+1,nstp)-v(i,j ,k,nstp)+ &
207 & v(i,j+1,k+1,nstp)-v(i,j+1,k,nstp))- &
208 & hz(i,j,k)*dv(i,k-1))
209 END DO
210 END DO
211 DO i=max(1,istr-1),min(iend+1,lm(ng))
212 dr(i,n(ng))=0.0_r8
213 du(i,n(ng))=0.0_r8
214 dv(i,n(ng))=0.0_r8
215 END DO
216 DO k=n(ng)-1,1,-1
217 DO i=max(1,istr-1),min(iend+1,lm(ng))
218 dr(i,k)=dr(i,k)-fc(i,k)*dr(i,k+1)
219 du(i,k)=du(i,k)-fc(i,k)*du(i,k+1)
220 dv(i,k)=dv(i,k)-fc(i,k)*dv(i,k+1)
221 END DO
222 END DO
223 DO k=1,n(ng)-1
224 DO i=max(1,istr-1),min(iend+1,lm(ng))
225 shear2=du(i,k)*du(i,k)+dv(i,k)*dv(i,k)
226 rig(i,j,k)=bvf(i,j,k)/(shear2+eps)
227!! Rig(i,j,k)=-gorho0*dR(i,k)/(shear2+eps)
228 END DO
229 END DO
230 END DO
231# else
232 DO k=1,n(ng)-1
233 DO j=max(1,jstr-1),min(jend+1,mm(ng))
234 DO i=max(1,istr-1),min(iend+1,lm(ng))
235 cff=0.5_r8/(z_r(i,j,k+1)-z_r(i,j,k))
236 shear2=(cff*(u(i ,j,k+1,nstp)-u(i ,j,k,nstp)+ &
237 & u(i+1,j,k+1,nstp)-u(i+1,j,k,nstp)))**2+ &
238 & (cff*(v(i,j ,k+1,nstp)-v(i,j ,k,nstp)+ &
239 & v(i,j+1,k+1,nstp)-v(i,j+1,k,nstp)))**2
240 rig(i,j,k)=bvf(i,j,k)/(shear2+eps)
241 END DO
242 END DO
243 END DO
244# endif
245# ifdef RI_HORAVG
246 DO k=1,n(ng)-1
247 IF (domain(ng)%Western_Edge(tile)) THEN
248 DO j=max(1,jstr-1),min(jend+1,mm(ng))
249 rig(istr-1,j,k)=rig(istr,j,k)
250 END DO
251 END IF
252 IF (domain(ng)%Eastern_Edge(tile)) THEN
253 DO j=max(1,jstr-1),min(jend+1,mm(ng))
254 rig(iend+1,j,k)=rig(iend,j,k)
255 END DO
256 END IF
257 IF (domain(ng)%Southern_Edge(tile)) THEN
258 DO i=max(1,istr-1),min(iend+1,lm(ng))
259 rig(i,jstr-1,k)=rig(i,jstr,k)
260 END DO
261 END IF
262 IF (domain(ng)%Northern_Edge(tile)) THEN
263 DO i=max(1,istr-1),min(iend+1,lm(ng))
264 rig(i,jend+1,k)=rig(i,jend,k)
265 END DO
266 END IF
267 IF (domain(ng)%SouthWest_Corner(tile)) THEN
268 rig(istr-1,jstr-1,k)=rig(istr,jstr,k)
269 END IF
270 IF (domain(ng)%NorthWest_Corner(tile)) THEN
271 rig(istr-1,jend+1,k)=rig(istr,jend,k)
272 END IF
273 IF (domain(ng)%SouthEast_Corner(tile)) THEN
274 rig(iend+1,jstr-1,k)=rig(iend,jstr,k)
275 END IF
276 IF (domain(ng)%NorthEast_Corner(tile)) THEN
277 rig(iend+1,jend+1,k)=rig(iend,jend,k)
278 END IF
279!
280! Smooth gradient Richardson number horizontally. Use Rig(:,:,0)
281! as scratch utility array.
282!
283 DO j=jstr-1,jend
284 DO i=istr-1,iend
285 rig(i,j,0)=0.25_r8*(rig(i,j ,k)+rig(i+1,j ,k)+ &
286 & rig(i,j+1,k)+rig(i+1,j+1,k))
287 END DO
288 END DO
289 DO j=jstr,jend
290 DO i=istr,iend
291 rig(i,j,k)=0.25_r8*(rig(i,j ,0)+rig(i-1,j ,0)+ &
292 & rig(i,j-1,0)+rig(i-1,j-1,0))
293 END DO
294 END DO
295 END DO
296# endif
297# ifdef RI_VERAVG
298!
299! Smooth gradient Richardson number vertically at the interior points.
300!
301 DO k=n(ng)-2,2,-1
302 DO j=jstr,jend
303 DO i=istr,iend
304 rig(i,j,k)=0.25_r8*rig(i,j,k-1)+ &
305 & 0.50_r8*rig(i,j,k )+ &
306 & 0.25_r8*rig(i,j,k+1)
307 END DO
308 END DO
309 END DO
310# endif
311# endif
312!
313!-----------------------------------------------------------------------
314! Compute "interior" viscosities and diffusivities everywhere as
315! the superposition of three processes: local Richardson number
316! instability due to resolved vertical shear, internal wave
317! breaking, and double diffusion.
318!-----------------------------------------------------------------------
319!
320 DO k=1,n(ng)-1
321 DO j=jstr,jend
322 DO i=istr,iend
323!
324! Compute interior diffusivity due to shear instability mixing.
325!
326# ifdef LMD_RIMIX
327 cff=min(1.0_r8,max(0.0_r8,rig(i,j,k))/lmd_ri0)
328 nu_sx=1.0_r8-cff*cff
329 nu_sx=nu_sx*nu_sx*nu_sx
330!
331! The shear mixing should be also a function of the actual magnitude
332! of the shear, see Polzin (1996, JPO, 1409-1425).
333!
334 shear2=bvf(i,j,k)/(rig(i,j,k)+eps)
335 cff=shear2*shear2/(shear2*shear2+16.0e-10_r8)
336 nu_sx=cff*nu_sx
337# else
338 nu_sx=0.0_r8
339# endif
340!
341! Compute interior diffusivity due to wave breaking (Gargett and
342! Holloway.
343!
344 cff=1.0_r8/sqrt(max(bvf(i,j,k),1.0e-7_r8))
345 lmd_iwm=1.0e-6_r8*cff
346 lmd_iws=1.0e-7_r8*cff
347! lmd_iwm=lmd_nuwm
348! lmd_iws=lmd_nuws
349!
350! Sum contributions due to internal wave breaking, shear instability
351! and convective diffusivity due to shear instability.
352!
353 akv(i,j,k)=lmd_iwm+lmd_nu0m*nu_sx
354 akt(i,j,k,itemp)=lmd_iws+lmd_nu0s*nu_sx
355# ifdef SALINITY
356 akt(i,j,k,isalt)=akt(i,j,k,itemp)
357# endif
358 END DO
359 END DO
360# ifdef LMD_DDMIX
361!
362!-----------------------------------------------------------------------
363! Compute double-diffusive mixing. It can occur when vertical
364! gradient of density is stable but the vertical gradient of
365! salinity (salt figering) or temperature (diffusive convection)
366! is unstable.
367!-----------------------------------------------------------------------
368!
369! Compute double-diffusive density ratio, Rrho.
370!
371 DO j=jstr,jend
372 DO i=istr,iend
373 dddt=t(i,j,k+1,nstp,itemp)-t(i,j,k,nstp,itemp)
374 ddds=t(i,j,k+1,nstp,isalt)-t(i,j,k,nstp,isalt)
375 ddds=sign(1.0_r8,ddds)*max(abs(ddds),1.0e-14_r8)
376 rrho=alfaobeta(i,j,k)*dddt/ddds
377!
378! Salt fingering case.
379!
380 IF ((rrho.gt.1.0_r8).and.(ddds.gt.0.0_r8)) THEN
381!
382! Compute interior diffusivity for double diffusive mixing of
383! salinity. Upper bound "Rrho" by "Rrho0"; (lmd_Rrho0=1.9,
384! lmd_nuf=0.001).
385!
386 rrho=min(rrho,lmd_rrho0)
387 nu_dds=1.0_r8-((rrho-1.0_r8)/(lmd_rrho0-1.0_r8))**2
388 nu_dds=lmd_nuf*nu_dds*nu_dds*nu_dds
389!
390! Compute interior diffusivity for double diffusive mixing
391! of temperature (lmd_fdd=0.7).
392!
393 nu_ddt=lmd_fdd*nu_dds
394!
395! Diffusive convection case.
396!
397 ELSE IF ((0.0_r8.lt.rrho).and.(rrho.lt.1.0_r8).and. &
398 & (ddds.lt.0.0_r8)) THEN
399!
400! Compute interior diffusivity for double diffusive mixing of
401! temperature (Marmorino and Caldwell, 1976); (lmd_nu=1.5e-6,
402! lmd_tdd1=0.909, lmd_tdd2=4.6, lmd_tdd3=0.54).
403!
404 nu_ddt=lmd_nu*lmd_tdd1* &
405 & exp(lmd_tdd2* &
406 & exp(-lmd_tdd3*((1.0_r8/rrho)-1.0_r8)))
407!
408! Compute interior diffusivity for double diffusive mixing
409! of salinity (lmd_sdd1=0.15, lmd_sdd2=1.85, lmd_sdd3=0.85).
410!
411 IF (rrho.lt.0.5_r8) THEN
412 nu_dds=nu_ddt*lmd_sdd1*rrho
413 ELSE
414 nu_dds=nu_ddt*(lmd_sdd2*rrho-lmd_sdd3)
415 END IF
416 ELSE
417 nu_ddt=0.0_r8
418 nu_dds=0.0_r8
419 END IF
420!
421! Add double diffusion contribution to temperature and salinity
422! mixing coefficients.
423!
424 akt(i,j,k,itemp)=akt(i,j,k,itemp)+nu_ddt
425# ifdef SALINITY
426 akt(i,j,k,isalt)=akt(i,j,k,isalt)+nu_dds
427# endif
428 END DO
429 END DO
430# endif
431 END DO
432
433 RETURN
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, dimension(:), allocatable mm
Definition mod_param.F:456
real(r8) lmd_tdd2
real(r8) lmd_tdd1
real(r8) lmd_sdd2
real(r8) lmd_ri0
real(r8) lmd_nu0s
real(r8) lmd_sdd1
real(r8) lmd_fdd
real(r8) lmd_rrho0
real(r8) lmd_nu
real(r8) lmd_sdd3
real(r8) lmd_nu0m
real(r8) lmd_nuf
real(r8) lmd_tdd3

References mod_param::domain, mod_scalars::isalt, mod_scalars::itemp, mod_param::lm, mod_scalars::lmd_fdd, mod_scalars::lmd_nu, mod_scalars::lmd_nu0m, mod_scalars::lmd_nu0s, mod_scalars::lmd_nuf, mod_scalars::lmd_ri0, mod_scalars::lmd_rrho0, mod_scalars::lmd_sdd1, mod_scalars::lmd_sdd2, mod_scalars::lmd_sdd3, mod_scalars::lmd_tdd1, mod_scalars::lmd_tdd2, mod_scalars::lmd_tdd3, and mod_param::mm.

Referenced by lmd_vmix().

Here is the caller graph for this function: