1 |
C $Header: /u/gcmpack/MITgcm/pkg/bling/bling_mixedlayer.F,v 1.1 2016/05/19 20:29:26 mmazloff Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "BLING_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
subroutine BLING_MIXEDLAYER( |
8 |
U sumMLDepth, |
9 |
I bi, bj, imin, imax, jmin, jmax, |
10 |
I myIter, myTime, myThid ) |
11 |
|
12 |
C ================================================================= |
13 |
C | subroutine bling_mixedlayer |
14 |
C | o Calculate mixed layer depth based on density criterion |
15 |
C ================================================================= |
16 |
|
17 |
implicit none |
18 |
|
19 |
C === Global variables === |
20 |
|
21 |
#include "SIZE.h" |
22 |
#include "EEPARAMS.h" |
23 |
#include "PARAMS.h" |
24 |
#include "FFIELDS.h" |
25 |
#include "GRID.h" |
26 |
#include "DYNVARS.h" |
27 |
#include "BLING_VARS.h" |
28 |
#include "PTRACERS_SIZE.h" |
29 |
#include "PTRACERS_PARAMS.h" |
30 |
#ifdef ALLOW_AUTODIFF |
31 |
# include "tamc.h" |
32 |
#endif |
33 |
|
34 |
C === Routine arguments === |
35 |
C bi,bj :: tile indices |
36 |
C iMin,iMax :: computation domain: 1rst index range |
37 |
C jMin,jMax :: computation domain: 2nd index range |
38 |
C myTime :: current time |
39 |
C myIter :: current timestep |
40 |
C myThid :: thread Id. number |
41 |
INTEGER bi, bj, imin, imax, jmin, jmax |
42 |
INTEGER myThid |
43 |
INTEGER myIter |
44 |
_RL myTime |
45 |
C === Output === |
46 |
C sumMLDepth :: mixed layer depth |
47 |
_RL sumMLDepth(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
48 |
|
49 |
C === Local variables === |
50 |
_RL dens_surf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
51 |
_RL dens_z (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
_RL delta_dens(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
53 |
integer i,j,k |
54 |
CEOP |
55 |
|
56 |
c --------------------------------------------------------------------- |
57 |
c Mixed layer depth |
58 |
|
59 |
DO j=jmin,jmax |
60 |
DO i=imin,imax |
61 |
SumMLDepth(i,j) = drf(1) |
62 |
ENDDO |
63 |
ENDDO |
64 |
|
65 |
c Surface density |
66 |
CALL FIND_RHO_2D( |
67 |
I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, |
68 |
I theta(1-OLx,1-OLy,1,bi,bj), salt(1-OLx,1-OLy,1,bi,bj), |
69 |
O dens_surf, |
70 |
I 1, bi, bj, myThid ) |
71 |
|
72 |
DO k=1,Nr |
73 |
DO j=jmin,jmax |
74 |
DO i=imin,imax |
75 |
if (k.eq.1) then |
76 |
delta_dens(i,j,1) = 0. _d 0 |
77 |
else |
78 |
delta_dens(i,j,k) = 9999. _d 0 |
79 |
endif |
80 |
ENDDO |
81 |
ENDDO |
82 |
ENDDO |
83 |
|
84 |
DO k = 2,Nr |
85 |
|
86 |
c Potential density |
87 |
CALL FIND_RHO_2D( |
88 |
I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, |
89 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
90 |
O dens_z, |
91 |
I k, bi, bj, myThid ) |
92 |
|
93 |
DO j=jmin,jmax |
94 |
DO i=imin,imax |
95 |
|
96 |
c SumMLDepth(i,j) = 0. _d 0 |
97 |
|
98 |
IF (hFacC(i,j,k,bi,bj) .gt. 0. _d 0) THEN |
99 |
delta_dens(i,j,k) = dens_z(i,j)-dens_surf(i,j) |
100 |
IF (delta_dens(i,j,k) .LT. 0.03 _d 0) THEN |
101 |
SumMLDepth(i,j) = SumMLDepth(i,j)+drF(k) |
102 |
ENDIF |
103 |
ENDIF |
104 |
|
105 |
ENDDO |
106 |
ENDDO |
107 |
ENDDO |
108 |
|
109 |
RETURN |
110 |
END |
111 |
|