1 |
jmc |
1.27 |
C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_implicit_r.F,v 1.26 2016/10/26 00:46:00 jmc Exp $ |
2 |
jmc |
1.24 |
C $Name: $ |
3 |
jmc |
1.1 |
|
4 |
|
|
#include "GAD_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
CBOP |
7 |
|
|
C !ROUTINE: GAD_IMPLICIT_R |
8 |
|
|
C !INTERFACE: |
9 |
jmc |
1.13 |
SUBROUTINE GAD_IMPLICIT_R( |
10 |
jmc |
1.1 |
I implicitAdvection, advectionScheme, tracerIdentity, |
11 |
jahn |
1.14 |
I deltaTLev, |
12 |
jmc |
1.19 |
I kappaRX, recip_hFac, wFld, tracer, |
13 |
jmc |
1.1 |
U gTracer, |
14 |
|
|
I bi, bj, myTime, myIter, myThid ) |
15 |
edhill |
1.4 |
C !DESCRIPTION: |
16 |
|
|
C Solve implicitly vertical advection and diffusion for one tracer. |
17 |
jmc |
1.1 |
|
18 |
|
|
C !USES: |
19 |
|
|
IMPLICIT NONE |
20 |
|
|
C == Global data == |
21 |
|
|
#include "SIZE.h" |
22 |
|
|
#include "EEPARAMS.h" |
23 |
|
|
#include "PARAMS.h" |
24 |
|
|
#include "GRID.h" |
25 |
jmc |
1.15 |
#include "SURFACE.h" |
26 |
jmc |
1.1 |
#include "GAD.h" |
27 |
|
|
|
28 |
edhill |
1.4 |
C !INPUT/OUTPUT PARAMETERS: |
29 |
|
|
C == Routine Arguments == |
30 |
|
|
C implicitAdvection :: if True, treat vertical advection implicitly |
31 |
|
|
C advectionScheme :: advection scheme to use |
32 |
|
|
C tracerIdentity :: Identifier for the tracer |
33 |
|
|
C kappaRX :: 3-D array for vertical diffusion coefficient |
34 |
jmc |
1.15 |
C recip_hFac :: inverse of cell open-depth factor |
35 |
jmc |
1.19 |
C wFld :: Advection velocity field, vertical component |
36 |
edhill |
1.4 |
C tracer :: tracer field at current time step |
37 |
|
|
C gTracer :: future tracer field |
38 |
|
|
C bi,bj :: tile indices |
39 |
|
|
C myTime :: current time |
40 |
|
|
C myIter :: current iteration number |
41 |
|
|
C myThid :: thread number |
42 |
jmc |
1.1 |
LOGICAL implicitAdvection |
43 |
|
|
INTEGER advectionScheme |
44 |
|
|
INTEGER tracerIdentity |
45 |
jahn |
1.14 |
_RL deltaTLev(Nr) |
46 |
jmc |
1.22 |
_RL kappaRX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
47 |
|
|
_RS recip_hFac(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
48 |
|
|
_RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
49 |
|
|
_RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
50 |
|
|
_RL gTracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
51 |
jmc |
1.1 |
INTEGER bi, bj |
52 |
|
|
_RL myTime |
53 |
|
|
INTEGER myIter, myThid |
54 |
|
|
|
55 |
jmc |
1.15 |
#ifdef ALLOW_DIAGNOSTICS |
56 |
|
|
C !FUNCTIONS: |
57 |
|
|
CHARACTER*4 GAD_DIAG_SUFX |
58 |
|
|
EXTERNAL GAD_DIAG_SUFX |
59 |
|
|
LOGICAL DIAGNOSTICS_IS_ON |
60 |
|
|
EXTERNAL DIAGNOSTICS_IS_ON |
61 |
|
|
#endif |
62 |
|
|
|
63 |
edhill |
1.4 |
C !LOCAL VARIABLES: |
64 |
|
|
C == Local variables == |
65 |
jmc |
1.15 |
C iMin,iMax,jMin,jMax :: computational domain |
66 |
|
|
C i,j,k :: loop indices |
67 |
|
|
C a5d :: 2nd lower diagonal of the pentadiagonal matrix |
68 |
|
|
C b5d :: 1rst lower diagonal of the pentadiagonal matrix |
69 |
|
|
C c5d :: main diagonal of the pentadiagonal matrix |
70 |
|
|
C d5d :: 1rst upper diagonal of the pentadiagonal matrix |
71 |
|
|
C e5d :: 2nd upper diagonal of the pentadiagonal matrix |
72 |
|
|
C rTrans :: vertical volume transport at interface k |
73 |
jmc |
1.24 |
C localTr :: local copy of tracer (for Non-Lin Adv.Scheme) |
74 |
edhill |
1.4 |
C diagonalNumber :: number of non-zero diagonals in the matrix |
75 |
|
|
C errCode :: > 0 if singular matrix |
76 |
jmc |
1.27 |
C msgBuf :: Informational/error message buffer |
77 |
jmc |
1.1 |
INTEGER iMin,iMax,jMin,jMax |
78 |
jmc |
1.15 |
PARAMETER( iMin = 1, iMax = sNx ) |
79 |
|
|
PARAMETER( jMin = 1, jMax = sNy ) |
80 |
jmc |
1.1 |
INTEGER i,j,k |
81 |
|
|
INTEGER diagonalNumber, errCode |
82 |
jmc |
1.22 |
_RL a5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
83 |
|
|
_RL b5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
84 |
|
|
_RL c5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
85 |
|
|
_RL d5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
86 |
|
|
_RL e5d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
87 |
jmc |
1.24 |
_RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
88 |
|
|
_RL localTr(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
89 |
jmc |
1.8 |
#ifdef ALLOW_DIAGNOSTICS |
90 |
|
|
CHARACTER*8 diagName |
91 |
jmc |
1.15 |
CHARACTER*4 diagSufx |
92 |
|
|
LOGICAL diagDif, diagAdv |
93 |
|
|
INTEGER km1, km2, kp1, kp2 |
94 |
jmc |
1.22 |
_RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
95 |
|
|
_RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
96 |
|
|
_RL div(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
97 |
|
|
_RL flx(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
98 |
jmc |
1.27 |
# ifdef SOLVE_DIAGONAL_LOWMEMORY |
99 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
100 |
|
|
# endif /* SOLVE_DIAGONAL_LOWMEMORY */ |
101 |
|
|
#endif /* ALLOW_DIAGNOSTICS */ |
102 |
jmc |
1.1 |
CEOP |
103 |
|
|
|
104 |
jmc |
1.7 |
C-- no need to solve anything with only 1 level: |
105 |
|
|
IF (Nr.GT.1) THEN |
106 |
jmc |
1.1 |
|
107 |
|
|
C-- Initialise |
108 |
|
|
DO k=1,Nr |
109 |
jmc |
1.22 |
DO j=1-OLy,sNy+OLy |
110 |
|
|
DO i=1-OLx,sNx+OLx |
111 |
jmc |
1.1 |
a5d(i,j,k) = 0. _d 0 |
112 |
|
|
b5d(i,j,k) = 0. _d 0 |
113 |
|
|
c5d(i,j,k) = 1. _d 0 |
114 |
|
|
d5d(i,j,k) = 0. _d 0 |
115 |
|
|
e5d(i,j,k) = 0. _d 0 |
116 |
|
|
ENDDO |
117 |
|
|
ENDDO |
118 |
|
|
ENDDO |
119 |
|
|
diagonalNumber = 1 |
120 |
|
|
|
121 |
|
|
IF (implicitDiffusion) THEN |
122 |
|
|
C-- set the tri-diagonal matrix to solve the implicit diffusion problem |
123 |
|
|
diagonalNumber = 3 |
124 |
|
|
C- 1rst lower diagonal : |
125 |
|
|
DO k=2,Nr |
126 |
|
|
DO j=jMin,jMax |
127 |
|
|
DO i=iMin,iMax |
128 |
jahn |
1.14 |
b5d(i,j,k) = -deltaTLev(k)*maskC(i,j,k-1,bi,bj) |
129 |
jmc |
1.15 |
& *recip_hFac(i,j,k)*recip_drF(k) |
130 |
jmc |
1.24 |
& *recip_deepFac2C(k)*recip_rhoFacC(k) |
131 |
jmc |
1.1 |
& *kappaRX(i,j, k )*recip_drC( k ) |
132 |
jmc |
1.24 |
& *deepFac2F( k )*rhoFacF( k ) |
133 |
jmc |
1.1 |
ENDDO |
134 |
|
|
ENDDO |
135 |
|
|
ENDDO |
136 |
|
|
C- 1rst upper diagonal : |
137 |
|
|
DO k=1,Nr-1 |
138 |
|
|
DO j=jMin,jMax |
139 |
|
|
DO i=iMin,iMax |
140 |
jahn |
1.14 |
d5d(i,j,k) = -deltaTLev(k)*maskC(i,j,k+1,bi,bj) |
141 |
jmc |
1.24 |
& *recip_hFac(i,j,k)*recip_drF(k) |
142 |
|
|
& *recip_deepFac2C(k)*recip_rhoFacC(k) |
143 |
|
|
& *kappaRX(i,j,k+1)*recip_drC(k+1) |
144 |
|
|
& *deepFac2F(k+1)*rhoFacF(k+1) |
145 |
jmc |
1.1 |
ENDDO |
146 |
|
|
ENDDO |
147 |
|
|
ENDDO |
148 |
|
|
C- Main diagonal : |
149 |
|
|
DO k=1,Nr |
150 |
|
|
DO j=jMin,jMax |
151 |
|
|
DO i=iMin,iMax |
152 |
jmc |
1.25 |
c5d(i,j,k) = 1. _d 0 - ( b5d(i,j,k) + d5d(i,j,k) ) |
153 |
|
|
C to recover older (prior to 2016-10-05) results: |
154 |
|
|
c c5d(i,j,k) = 1. _d 0 - b5d(i,j,k) - d5d(i,j,k) |
155 |
jmc |
1.1 |
ENDDO |
156 |
|
|
ENDDO |
157 |
|
|
ENDDO |
158 |
|
|
|
159 |
|
|
C-- end if implicitDiffusion |
160 |
|
|
ENDIF |
161 |
|
|
|
162 |
|
|
IF (implicitAdvection) THEN |
163 |
|
|
|
164 |
jmc |
1.24 |
C-- Non-Linear Advection scheme: keep a local copy of tracer field |
165 |
|
|
IF ( advectionScheme.EQ.ENUM_FLUX_LIMIT .OR. |
166 |
|
|
& advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN |
167 |
|
|
IF ( multiDimAdvection ) THEN |
168 |
|
|
DO k=1,Nr |
169 |
|
|
DO j=1-OLy,sNy+OLy |
170 |
|
|
DO i=1-OLx,sNx+OLx |
171 |
|
|
localTr(i,j,k) = gTracer(i,j,k) |
172 |
|
|
ENDDO |
173 |
|
|
ENDDO |
174 |
|
|
ENDDO |
175 |
|
|
ELSE |
176 |
|
|
DO k=1,Nr |
177 |
|
|
DO j=1-OLy,sNy+OLy |
178 |
|
|
DO i=1-OLx,sNx+OLx |
179 |
|
|
localTr(i,j,k) = tracer(i,j,k) |
180 |
|
|
ENDDO |
181 |
|
|
ENDDO |
182 |
|
|
ENDDO |
183 |
|
|
ENDIF |
184 |
|
|
ENDIF |
185 |
|
|
|
186 |
jmc |
1.2 |
DO k=Nr,1,-1 |
187 |
|
|
|
188 |
|
|
C-- Compute transport |
189 |
|
|
IF (k.EQ.1) THEN |
190 |
jmc |
1.22 |
DO j=1-OLy,sNy+OLy |
191 |
|
|
DO i=1-OLx,sNx+OLx |
192 |
jmc |
1.13 |
rTrans(i,j) = 0. _d 0 |
193 |
jmc |
1.2 |
ENDDO |
194 |
|
|
ENDDO |
195 |
|
|
ELSE |
196 |
jmc |
1.22 |
DO j=1-OLy,sNy+OLy |
197 |
|
|
DO i=1-OLx,sNx+OLx |
198 |
jmc |
1.19 |
rTrans(i,j) = wFld(i,j,k)*rA(i,j,bi,bj) |
199 |
jmc |
1.24 |
& *deepFac2F(k)*rhoFacF(k) |
200 |
jmc |
1.19 |
& *maskC(i,j,k-1,bi,bj) |
201 |
jmc |
1.1 |
ENDDO |
202 |
jmc |
1.2 |
ENDDO |
203 |
|
|
ENDIF |
204 |
|
|
|
205 |
jmc |
1.5 |
#ifdef ALLOW_AIM |
206 |
|
|
C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr |
207 |
jmc |
1.15 |
IF ( k.GE.2 .AND. |
208 |
|
|
& (.NOT.useAIM .OR.tracerIdentity.NE.GAD_SALINITY .OR.k.LT.Nr) |
209 |
jmc |
1.5 |
& ) THEN |
210 |
|
|
#else |
211 |
jmc |
1.15 |
IF ( k.GE.2 ) THEN |
212 |
jmc |
1.5 |
#endif |
213 |
jmc |
1.1 |
|
214 |
jmc |
1.2 |
IF ( advectionScheme.EQ.ENUM_CENTERED_2ND ) THEN |
215 |
|
|
diagonalNumber = 3 |
216 |
|
|
CALL GAD_C2_IMPL_R( bi,bj,k, iMin,iMax,jMin,jMax, |
217 |
jmc |
1.15 |
I deltaTLev, rTrans, recip_hFac, |
218 |
jmc |
1.2 |
U b5d, c5d, d5d, |
219 |
jmc |
1.11 |
I myThid ) |
220 |
|
|
ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST |
221 |
|
|
& .OR. advectionScheme.EQ.ENUM_DST2 ) THEN |
222 |
|
|
diagonalNumber = 3 |
223 |
|
|
CALL GAD_DST2U1_IMPL_R( bi,bj,k, iMin,iMax,jMin,jMax, |
224 |
jmc |
1.15 |
I advectionScheme, deltaTLev, |
225 |
|
|
I rTrans, recip_hFac, |
226 |
jmc |
1.11 |
U b5d, c5d, d5d, |
227 |
|
|
I myThid ) |
228 |
|
|
ELSEIF ( advectionScheme.EQ.ENUM_FLUX_LIMIT ) THEN |
229 |
jmc |
1.2 |
diagonalNumber = 3 |
230 |
|
|
CALL GAD_FLUXLIMIT_IMPL_R( bi,bj,k, iMin,iMax,jMin,jMax, |
231 |
jmc |
1.24 |
I deltaTLev, rTrans, recip_hFac, localTr, |
232 |
jmc |
1.2 |
U b5d, c5d, d5d, |
233 |
jmc |
1.11 |
I myThid ) |
234 |
|
|
ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_3RD |
235 |
|
|
& .OR. advectionScheme.EQ.ENUM_CENTERED_4TH |
236 |
|
|
& .OR. advectionScheme.EQ.ENUM_DST3 ) THEN |
237 |
jmc |
1.2 |
diagonalNumber = 5 |
238 |
|
|
CALL GAD_U3C4_IMPL_R( bi,bj,k, iMin,iMax,jMin,jMax, |
239 |
jmc |
1.15 |
I advectionScheme, deltaTLev, |
240 |
|
|
I rTrans, recip_hFac, |
241 |
jmc |
1.2 |
U a5d, b5d, c5d, d5d, e5d, |
242 |
jmc |
1.11 |
I myThid ) |
243 |
|
|
ELSEIF ( advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN |
244 |
|
|
diagonalNumber = 5 |
245 |
|
|
CALL GAD_DST3FL_IMPL_R( bi,bj,k, iMin,iMax,jMin,jMax, |
246 |
jmc |
1.24 |
I deltaTLev, rTrans, recip_hFac, localTr, |
247 |
jmc |
1.11 |
U a5d, b5d, c5d, d5d, e5d, |
248 |
|
|
I myThid ) |
249 |
jmc |
1.2 |
ELSE |
250 |
|
|
STOP 'GAD_IMPLICIT_R: Adv.Scheme in Impl form not yet coded' |
251 |
|
|
ENDIF |
252 |
|
|
|
253 |
|
|
ENDIF |
254 |
jmc |
1.1 |
|
255 |
jmc |
1.2 |
C-- end k loop |
256 |
|
|
ENDDO |
257 |
jmc |
1.1 |
|
258 |
|
|
C-- end if implicitAdvection |
259 |
|
|
ENDIF |
260 |
|
|
|
261 |
|
|
IF ( diagonalNumber .EQ. 3 ) THEN |
262 |
|
|
C-- Solve tri-diagonal system : |
263 |
jmc |
1.26 |
errCode = -1 |
264 |
jmc |
1.1 |
CALL SOLVE_TRIDIAGONAL( iMin,iMax, jMin,jMax, |
265 |
|
|
I b5d, c5d, d5d, |
266 |
|
|
U gTracer, |
267 |
|
|
O errCode, |
268 |
|
|
I bi, bj, myThid ) |
269 |
|
|
IF (errCode.GE.1) THEN |
270 |
|
|
STOP 'GAD_IMPLICIT_R: error when solving 3-Diag problem' |
271 |
jmc |
1.2 |
ENDIF |
272 |
|
|
ELSEIF ( diagonalNumber .EQ. 5 ) THEN |
273 |
|
|
C-- Solve penta-diagonal system : |
274 |
jmc |
1.26 |
errCode = -1 |
275 |
jmc |
1.2 |
CALL SOLVE_PENTADIAGONAL( iMin,iMax, jMin,jMax, |
276 |
|
|
I a5d, b5d, c5d, d5d, e5d, |
277 |
|
|
U gTracer, |
278 |
|
|
O errCode, |
279 |
|
|
I bi, bj, myThid ) |
280 |
|
|
IF (errCode.GE.1) THEN |
281 |
|
|
STOP 'GAD_IMPLICIT_R: error when solving 5-Diag problem' |
282 |
jmc |
1.1 |
ENDIF |
283 |
|
|
ELSEIF ( diagonalNumber .NE. 1 ) THEN |
284 |
|
|
STOP 'GAD_IMPLICIT_R: no solver available' |
285 |
|
|
ENDIF |
286 |
|
|
|
287 |
jmc |
1.8 |
#ifdef ALLOW_DIAGNOSTICS |
288 |
|
|
C-- Set diagnostic suffix for the current tracer |
289 |
jmc |
1.15 |
IF ( useDiagnostics ) THEN |
290 |
jmc |
1.8 |
diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid ) |
291 |
|
|
diagName = 'DFrI'//diagSufx |
292 |
jmc |
1.15 |
diagDif = implicitDiffusion |
293 |
|
|
IF ( diagDif ) diagDif = DIAGNOSTICS_IS_ON(diagName,myThid) |
294 |
|
|
diagName = 'ADVr'//diagSufx |
295 |
|
|
diagAdv = implicitAdvection |
296 |
|
|
IF ( diagAdv ) diagAdv = DIAGNOSTICS_IS_ON(diagName,myThid) |
297 |
|
|
|
298 |
|
|
IF ( diagDif .OR. diagAdv ) THEN |
299 |
|
|
DO j=1-OLy,sNy+OLy |
300 |
|
|
DO i=1-OLx,sNx+OLx |
301 |
gforget |
1.17 |
flx(i,j) = 0. _d 0 |
302 |
jmc |
1.15 |
ENDDO |
303 |
|
|
ENDDO |
304 |
jmc |
1.27 |
C-- start diagnostics k loop |
305 |
jmc |
1.15 |
DO k= Nr,1,-1 |
306 |
jmc |
1.27 |
|
307 |
jmc |
1.15 |
IF ( implicitDiffusion .AND. k.GE.2 ) THEN |
308 |
|
|
DO j=jMin,jMax |
309 |
|
|
DO i=iMin,iMax |
310 |
|
|
df(i,j) = |
311 |
heimbach |
1.18 |
cc#ifdef ALLOW_AUTODIFF_OPENAD |
312 |
|
|
cc & -rA(i,j,bi,bj)%v |
313 |
|
|
cc#else |
314 |
jmc |
1.24 |
& -rA(i,j,bi,bj)*deepFac2F(k)*rhoFacF(k) |
315 |
heimbach |
1.18 |
cc#endif |
316 |
jmc |
1.24 |
& * kappaRX(i,j,k)*recip_drC(k)*rkSign |
317 |
jmc |
1.22 |
& * (gTracer(i,j,k) - gTracer(i,j,k-1)) |
318 |
jmc |
1.15 |
& * maskC(i,j,k,bi,bj) |
319 |
|
|
& * maskC(i,j,k-1,bi,bj) |
320 |
|
|
ENDDO |
321 |
|
|
ENDDO |
322 |
|
|
ELSE |
323 |
jmc |
1.8 |
DO j=1-OLy,sNy+OLy |
324 |
|
|
DO i=1-OLx,sNx+OLx |
325 |
|
|
df(i,j) = 0. _d 0 |
326 |
|
|
ENDDO |
327 |
|
|
ENDDO |
328 |
jmc |
1.15 |
ENDIF |
329 |
jmc |
1.27 |
|
330 |
jmc |
1.15 |
C- Note: Needs to explicitly increment counter (call DIAGNOSTICS_COUNT) |
331 |
|
|
C since skipping k=1 DIAGNOSTICS_FILL call. |
332 |
|
|
IF ( diagDif .AND. k.GE.2 ) THEN |
333 |
|
|
diagName = 'DFrI'//diagSufx |
334 |
|
|
CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid) |
335 |
|
|
IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid) |
336 |
rpa |
1.20 |
#ifdef ALLOW_LAYERS |
337 |
rpa |
1.21 |
IF ( useLayers ) THEN |
338 |
rpa |
1.23 |
CALL LAYERS_FILL( df, tracerIdentity, 'DFR', |
339 |
rpa |
1.20 |
& k, 1, 2,bi,bj, myThid ) |
340 |
jmc |
1.22 |
ENDIF |
341 |
rpa |
1.20 |
#endif /* ALLOW_LAYERS */ |
342 |
jmc |
1.15 |
ENDIF |
343 |
jmc |
1.27 |
|
344 |
jmc |
1.15 |
IF ( diagAdv ) THEN |
345 |
jmc |
1.27 |
#ifdef SOLVE_DIAGONAL_LOWMEMORY |
346 |
|
|
diagName = 'ADVr'//diagSufx |
347 |
|
|
WRITE(msgBuf,'(4A)') 'GAD_IMPLICIT_R: ', |
348 |
|
|
& 'unable to compute Diagnostic "', diagName, '" with' |
349 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
350 |
|
|
& SQUEEZE_RIGHT, myThid ) |
351 |
|
|
WRITE(msgBuf,'(4A)') 'GAD_IMPLICIT_R: ', |
352 |
|
|
& '#define SOLVE_DIAGONAL_LOWMEMORY (in CPP_OPTIONS.h)' |
353 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
354 |
|
|
& SQUEEZE_RIGHT, myThid ) |
355 |
|
|
STOP 'ABNORMAL END: S/R GAD_IMPLICIT_R' |
356 |
|
|
#endif /* SOLVE_DIAGONAL_LOWMEMORY */ |
357 |
jmc |
1.15 |
km1=MAX(1,k-1) |
358 |
|
|
km2=MAX(1,k-2) |
359 |
|
|
kp1=MIN(Nr,k+1) |
360 |
|
|
kp2=MIN(Nr,k+2) |
361 |
|
|
C-- Flux_divergence*deltaT = Tr^n - Tr^n+1 = [A-I](Tr^n+1) |
362 |
|
|
C = deltaT*rkSign*[ Flx_k+1 - Flx_k ]/dz |
363 |
|
|
DO j=jMin,jMax |
364 |
|
|
DO i=iMin,iMax |
365 |
jmc |
1.22 |
div(i,j) = gTracer(i,j,k)*( c5d(i,j,k) - 1. _d 0 ) |
366 |
|
|
& + gTracer(i,j,km1)*b5d(i,j,k) |
367 |
|
|
& + gTracer(i,j,kp1)*d5d(i,j,k) |
368 |
jmc |
1.15 |
ENDDO |
369 |
|
|
ENDDO |
370 |
|
|
IF ( diagonalNumber .EQ. 5 ) THEN |
371 |
|
|
DO j=jMin,jMax |
372 |
|
|
DO i=iMin,iMax |
373 |
|
|
div(i,j) = div(i,j) |
374 |
jmc |
1.22 |
& + gTracer(i,j,km2)*a5d(i,j,k) |
375 |
|
|
& + gTracer(i,j,kp2)*e5d(i,j,k) |
376 |
jmc |
1.15 |
ENDDO |
377 |
|
|
ENDDO |
378 |
|
|
ENDIF |
379 |
|
|
#ifdef NONLIN_FRSURF |
380 |
|
|
IF ( nonlinFreeSurf.GT.0 ) THEN |
381 |
|
|
C-- use future hFac to stay consistent with solver matrix |
382 |
|
|
IF ( select_rStar.GT.0 ) THEN |
383 |
|
|
DO j=jMin,jMax |
384 |
|
|
DO i=iMin,iMax |
385 |
|
|
div(i,j) = div(i,j)*h0FacC(i,j,k,bi,bj)*drF(k) |
386 |
|
|
& *rStarFacC(i,j,bi,bj) |
387 |
|
|
ENDDO |
388 |
|
|
ENDDO |
389 |
|
|
ELSEIF ( selectSigmaCoord.NE.0 ) THEN |
390 |
|
|
DO j=jMin,jMax |
391 |
|
|
DO i=iMin,iMax |
392 |
|
|
div(i,j) = div(i,j)*( |
393 |
|
|
& _hFacC(i,j,k,bi,bj)*drF(k) |
394 |
jmc |
1.22 |
& + dBHybSigF(k)*dEtaHdt(i,j,bi,bj)*deltaTFreeSurf |
395 |
jmc |
1.15 |
& ) |
396 |
|
|
ENDDO |
397 |
|
|
ENDDO |
398 |
|
|
ELSE |
399 |
|
|
DO j=jMin,jMax |
400 |
|
|
DO i=iMin,iMax |
401 |
|
|
IF ( k.EQ.kSurfC(i,j,bi,bj) ) THEN |
402 |
|
|
div(i,j) = div(i,j)*hFac_surfC(i,j,bi,bj)*drF(k) |
403 |
|
|
ELSE |
404 |
|
|
div(i,j) = div(i,j)*_hFacC(i,j,k,bi,bj)*drF(k) |
405 |
|
|
ENDIF |
406 |
|
|
ENDDO |
407 |
|
|
ENDDO |
408 |
|
|
ENDIF |
409 |
|
|
ELSE |
410 |
|
|
#else /* NONLIN_FRSURF */ |
411 |
|
|
IF ( .TRUE. ) THEN |
412 |
|
|
#endif /* NONLIN_FRSURF */ |
413 |
|
|
C-- use current hFac (consistent with solver matrix) |
414 |
|
|
DO j=jMin,jMax |
415 |
|
|
DO i=iMin,iMax |
416 |
|
|
div(i,j) = div(i,j)*_hFacC(i,j,k,bi,bj)*drF(k) |
417 |
|
|
ENDDO |
418 |
jmc |
1.8 |
ENDDO |
419 |
jmc |
1.15 |
ENDIF |
420 |
|
|
DO j=jMin,jMax |
421 |
|
|
DO i=iMin,iMax |
422 |
gforget |
1.17 |
flx(i,j) = flx(i,j) |
423 |
heimbach |
1.18 |
cc#ifdef ALLOW_AUTODIFF_OPENAD |
424 |
|
|
cc & - rkSign*div(i,j)*rA(i,j,bi,bj)%v/deltaTLev(k) |
425 |
|
|
cc#else |
426 |
jmc |
1.24 |
& - rkSign*div(i,j)*rA(i,j,bi,bj) |
427 |
|
|
& *deepFac2C(k)*rhoFacC(k)/deltaTLev(k) |
428 |
heimbach |
1.18 |
cc#endif |
429 |
gforget |
1.17 |
af(i,j) = flx(i,j) - df(i,j) |
430 |
jmc |
1.8 |
ENDDO |
431 |
jmc |
1.15 |
ENDDO |
432 |
|
|
diagName = 'ADVr'//diagSufx |
433 |
|
|
CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid) |
434 |
rpa |
1.23 |
#ifdef ALLOW_LAYERS |
435 |
|
|
IF ( useLayers ) THEN |
436 |
|
|
CALL LAYERS_FILL(af,tracerIdentity,'AFR', |
437 |
|
|
& k,1,2,bi,bj,myThid) |
438 |
|
|
ENDIF |
439 |
|
|
#endif /* ALLOW_LAYERS */ |
440 |
jmc |
1.8 |
ENDIF |
441 |
jmc |
1.27 |
|
442 |
|
|
C-- end diagnostics k loop |
443 |
jmc |
1.8 |
ENDDO |
444 |
|
|
ENDIF |
445 |
|
|
ENDIF |
446 |
|
|
#endif /* ALLOW_DIAGNOSTICS */ |
447 |
|
|
|
448 |
jmc |
1.7 |
C-- end if Nr > 1 |
449 |
|
|
ENDIF |
450 |
|
|
|
451 |
jmc |
1.1 |
RETURN |
452 |
|
|
END |