1 |
C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_rtransport.F,v 1.12 2003/01/21 19:34:13 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "GMREDI_OPTIONS.h" |
5 |
|
6 |
subroutine GMREDI_RTRANSPORT( |
7 |
I iMin,iMax,jMin,jMax,bi,bj,K, |
8 |
I Tracer,tracerIdentity, |
9 |
U df, |
10 |
I myThid) |
11 |
C /==========================================================\ |
12 |
C | o SUBROUTINE GMREDI_RTRANSPORT | |
13 |
C | Add vertical transport terms from GM/Redi | |
14 |
C | parameterization. | |
15 |
C |==========================================================| |
16 |
C \==========================================================/ |
17 |
IMPLICIT NONE |
18 |
|
19 |
C == GLobal variables == |
20 |
#include "SIZE.h" |
21 |
#include "EEPARAMS.h" |
22 |
#include "PARAMS.h" |
23 |
#include "GRID.h" |
24 |
#include "GMREDI.h" |
25 |
|
26 |
#ifdef ALLOW_AUTODIFF_TAMC |
27 |
#include "tamc.h" |
28 |
#include "tamc_keys.h" |
29 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
30 |
|
31 |
C == Routine arguments == |
32 |
C iMin,iMax,jMin, - Range of points for which calculation |
33 |
C jMax,bi,bj,k results will be set. |
34 |
C xA - Area of X face |
35 |
C Tracer - 3D Tracer field |
36 |
C df - Diffusive flux component work array. |
37 |
INTEGER iMin,iMax,jMin,jMax,bi,bj,k |
38 |
_RL Tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
39 |
INTEGER tracerIdentity |
40 |
_RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
41 |
INTEGER myThid |
42 |
|
43 |
#ifdef ALLOW_GMREDI |
44 |
|
45 |
C == Local variables == |
46 |
C I, J - Loop counters |
47 |
INTEGER I, J |
48 |
_RL dTdx (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
49 |
_RL dTdy (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
50 |
#ifdef GM_BOLUS_ADVEC |
51 |
_RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
#endif |
53 |
|
54 |
#ifdef ALLOW_AUTODIFF_TAMC |
55 |
act0 = tracerIdentity - 1 |
56 |
max0 = maxpass |
57 |
act1 = bi - myBxLo(myThid) |
58 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
59 |
act2 = bj - myByLo(myThid) |
60 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
61 |
act3 = myThid - 1 |
62 |
max3 = nTx*nTy |
63 |
act4 = ikey_dynamics - 1 |
64 |
igadkey = (act0 + 1) |
65 |
& + act1*max0 |
66 |
& + act2*max0*max1 |
67 |
& + act3*max0*max1*max2 |
68 |
& + act4*max0*max1*max2*max3 |
69 |
kkey = (igadkey-1)*Nr + k |
70 |
if (tracerIdentity.GT.maxpass) then |
71 |
print *, 'ph-pass gmredi_rtrans ', maxpass, tracerIdentity |
72 |
STOP 'maxpass seems smaller than tracerIdentity' |
73 |
endif |
74 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
75 |
|
76 |
C Surface flux is zero |
77 |
IF (useGMRedi .AND. K.GT.1) THEN |
78 |
|
79 |
C- Horizontal gradients interpolated to W points |
80 |
DO j=jMin,jMax |
81 |
DO i=iMin,iMax |
82 |
dTdx(i,j) = op5*( |
83 |
& +op5*(_maskW(i+1,j,k,bi,bj) |
84 |
& *_recip_dxC(i+1,j,bi,bj)* |
85 |
& (Tracer(i+1,j,k,bi,bj)-Tracer(i,j,k,bi,bj)) |
86 |
& +_maskW(i,j,k,bi,bj) |
87 |
& *_recip_dxC(i,j,bi,bj)* |
88 |
& (Tracer(i,j,k,bi,bj)-Tracer(i-1,j,k,bi,bj))) |
89 |
& +op5*(_maskW(i+1,j,k-1,bi,bj) |
90 |
& *_recip_dxC(i+1,j,bi,bj)* |
91 |
& (Tracer(i+1,j,k-1,bi,bj)-Tracer(i,j,k-1,bi,bj)) |
92 |
& +_maskW(i,j,k-1,bi,bj) |
93 |
& *_recip_dxC(i,j,bi,bj)* |
94 |
& (Tracer(i,j,k-1,bi,bj)-Tracer(i-1,j,k-1,bi,bj))) |
95 |
& ) |
96 |
|
97 |
dTdy(i,j) = op5*( |
98 |
& +op5*(_maskS(i,j,k,bi,bj) |
99 |
& *_recip_dyC(i,j,bi,bj)* |
100 |
& (Tracer(i,j,k,bi,bj)-Tracer(i,j-1,k,bi,bj)) |
101 |
& +_maskS(i,j+1,k,bi,bj) |
102 |
& *_recip_dyC(i,j+1,bi,bj)* |
103 |
& (Tracer(i,j+1,k,bi,bj)-Tracer(i,j,k,bi,bj))) |
104 |
& +op5*(_maskS(i,j,k-1,bi,bj) |
105 |
& *_recip_dyC(i,j,bi,bj)* |
106 |
& (Tracer(i,j,k-1,bi,bj)-Tracer(i,j-1,k-1,bi,bj)) |
107 |
& +_maskS(i,j+1,k-1,bi,bj) |
108 |
& *_recip_dyC(i,j+1,bi,bj)* |
109 |
& (Tracer(i,j+1,k-1,bi,bj)-Tracer(i,j,k-1,bi,bj))) |
110 |
& ) |
111 |
ENDDO |
112 |
ENDDO |
113 |
|
114 |
#ifdef GM_AUTODIFF_EXCESSIVE_STORE |
115 |
CADJ STORE dTdx(:,:) = |
116 |
CADJ & comlev1_gmredi_k_gad, key=kkey, byte=isbyte |
117 |
CADJ STORE dTdy(:,:) = |
118 |
CADJ & comlev1_gmredi_k_gad, key=kkey, byte=isbyte |
119 |
#endif |
120 |
|
121 |
C- Off-diagonal components of vertical flux |
122 |
DO j=jMin,jMax |
123 |
DO i=iMin,iMax |
124 |
df(i,j) = df(i,j) |
125 |
& - _rA(i,j,bi,bj) |
126 |
& *( Kwx(i,j,k,bi,bj)*dTdx(i,j)+Kwy(i,j,k,bi,bj)*dTdy(i,j) ) |
127 |
|
128 |
ENDDO |
129 |
ENDDO |
130 |
|
131 |
#ifdef GM_BOLUS_ADVEC |
132 |
IF (GM_AdvForm .AND. GM_AdvSeparate) THEN |
133 |
DO j=jMin,jMax |
134 |
DO i=iMin,iMax |
135 |
rTrans(i,j) = |
136 |
& dyG(i+1,j,bi,bj)*GM_PsiX(i+1,j,k,bi,bj) |
137 |
& -dyG( i ,j,bi,bj)*GM_PsiX( i ,j,k,bi,bj) |
138 |
& +dxG(i,j+1,bi,bj)*GM_PsiY(i,j+1,k,bi,bj) |
139 |
& -dxG(i, j ,bi,bj)*GM_PsiY(i, j ,k,bi,bj) |
140 |
ENDDO |
141 |
ENDDO |
142 |
#ifdef GM_AUTODIFF_EXCESSIVE_STORE |
143 |
CADJ STORE rtrans(:,:) = |
144 |
CADJ & comlev1_gmredi_k_gad, key=kkey, byte=isbyte |
145 |
#endif |
146 |
DO j=jMin,jMax |
147 |
DO i=iMin,iMax |
148 |
df(i,j) = df(i,j) |
149 |
& +rTrans(i,j)*op5 |
150 |
& *(Tracer(i,j,k,bi,bj)+Tracer(i,j,k-1,bi,bj)) |
151 |
ENDDO |
152 |
ENDDO |
153 |
ENDIF |
154 |
#endif /* GM_BOLUS_ADVEC */ |
155 |
|
156 |
c IF (.NOT.implicitDiffusion) THEN |
157 |
c |
158 |
c This vertical diffusion term is currently implemented |
159 |
c by adding the VisbeckK*Kwz diffusivity to KappaRT/S |
160 |
c See calc_diffusivity.F and calc_gt.F (calc_gs.F) |
161 |
c |
162 |
c DO j=jMin,jMax |
163 |
c DO i=iMin,iMax |
164 |
c df(i,j) = df(i,j) - _rA(i,j,bi,bj) |
165 |
c & *maskUp(i,j)*VisbeckK(i,j,bi,bj)*Kwz(i,j,k,bi,bj) |
166 |
c & *recip_drC(k)*rkfac |
167 |
c & *(Tracer(i,j,k-1,bi,bj)-Tracer(i,j,k,bi,bj)) |
168 |
c ENDDO |
169 |
c ENDDO |
170 |
c ENDIF |
171 |
|
172 |
ENDIF |
173 |
#endif /* ALLOW_GMREDI */ |
174 |
|
175 |
RETURN |
176 |
END |