/[MITgcm]/MITgcm/model/src/do_oceanic_phys.F
ViewVC logotype

Contents of /MITgcm/model/src/do_oceanic_phys.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.126 - (show annotations) (download)
Sat Apr 13 20:47:18 2013 UTC (11 years, 1 month ago) by heimbach
Branch: MAIN
Changes since 1.125: +9 -3 lines
Add headers for THSICE adjoint.

1 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.125 2013/03/30 01:25:44 heimbach Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_AUTODIFF_TAMC
8 # ifdef ALLOW_GMREDI
9 # include "GMREDI_OPTIONS.h"
10 # endif
11 # ifdef ALLOW_KPP
12 # include "KPP_OPTIONS.h"
13 # endif
14 # ifdef ALLOW_SEAICE
15 # include "SEAICE_OPTIONS.h"
16 # endif
17 # ifdef ALLOW_EXF
18 # include "EXF_OPTIONS.h"
19 # endif
20 #endif /* ALLOW_AUTODIFF_TAMC */
21
22 CBOP
23 C !ROUTINE: DO_OCEANIC_PHYS
24 C !INTERFACE:
25 SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)
26 C !DESCRIPTION: \bv
27 C *==========================================================*
28 C | SUBROUTINE DO_OCEANIC_PHYS
29 C | o Controlling routine for oceanic physics and
30 C | parameterization
31 C *==========================================================*
32 C | o originally, part of S/R thermodynamics
33 C *==========================================================*
34 C \ev
35
36 C !USES:
37 IMPLICIT NONE
38 C == Global variables ===
39 #include "SIZE.h"
40 #include "EEPARAMS.h"
41 #include "PARAMS.h"
42 #include "GRID.h"
43 #include "DYNVARS.h"
44 #ifdef ALLOW_TIMEAVE
45 #include "TIMEAVE_STATV.h"
46 #endif
47 #ifdef ALLOW_BALANCE_FLUXES
48 #include "FFIELDS.h"
49 #endif
50
51 #ifdef ALLOW_AUTODIFF_TAMC
52 # include "AUTODIFF_MYFIELDS.h"
53 # include "tamc.h"
54 # include "tamc_keys.h"
55 #ifndef ALLOW_BALANCE_FLUXES
56 # include "FFIELDS.h"
57 #endif
58 # include "SURFACE.h"
59 # include "EOS.h"
60 # ifdef ALLOW_KPP
61 # include "KPP.h"
62 # endif
63 # ifdef ALLOW_GGL90
64 # include "GGL90.h"
65 # endif
66 # ifdef ALLOW_GMREDI
67 # include "GMREDI.h"
68 # endif
69 # ifdef ALLOW_EBM
70 # include "EBM.h"
71 # endif
72 # ifdef ALLOW_EXF
73 # include "ctrl.h"
74 # include "EXF_FIELDS.h"
75 # ifdef ALLOW_BULKFORMULAE
76 # include "EXF_CONSTANTS.h"
77 # endif
78 # endif
79 # ifdef ALLOW_SEAICE
80 # include "SEAICE_SIZE.h"
81 # include "SEAICE.h"
82 # include "SEAICE_PARAMS.h"
83 # endif
84 # ifdef ALLOW_THSICE
85 # include "THSICE_VARS.h"
86 # endif
87 # ifdef ALLOW_SALT_PLUME
88 # include "SALT_PLUME.h"
89 # endif
90 #endif /* ALLOW_AUTODIFF_TAMC */
91
92 C !INPUT/OUTPUT PARAMETERS:
93 C == Routine arguments ==
94 C myTime :: Current time in simulation
95 C myIter :: Current iteration number in simulation
96 C myThid :: Thread number for this instance of the routine.
97 _RL myTime
98 INTEGER myIter
99 INTEGER myThid
100
101 C !LOCAL VARIABLES:
102 C == Local variables
103 C rhoK, rhoKm1 :: Density at current level, and level above
104 C iMin, iMax :: Ranges and sub-block indices on which calculations
105 C jMin, jMax are applied.
106 C bi, bj :: tile indices
107 C i,j,k :: loop indices
108 _RL rhoKp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109 _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
111 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
112 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
113 INTEGER iMin, iMax
114 INTEGER jMin, jMax
115 INTEGER bi, bj
116 INTEGER i, j, k
117 INTEGER doDiagsRho
118 #ifdef ALLOW_DIAGNOSTICS
119 LOGICAL DIAGNOSTICS_IS_ON
120 EXTERNAL DIAGNOSTICS_IS_ON
121 #endif /* ALLOW_DIAGNOSTICS */
122
123 CEOP
124
125 #ifdef ALLOW_AUTODIFF_TAMC
126 C-- dummy statement to end declaration part
127 itdkey = 1
128 #endif /* ALLOW_AUTODIFF_TAMC */
129
130 #ifdef ALLOW_DEBUG
131 IF (debugMode) CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
132 #endif
133
134 doDiagsRho = 0
135 #ifdef ALLOW_DIAGNOSTICS
136 IF ( useDiagnostics .AND. fluidIsWater ) THEN
137 IF ( DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) )
138 & doDiagsRho = doDiagsRho + 1
139 IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) )
140 & doDiagsRho = doDiagsRho + 2
141 IF ( DIAGNOSTICS_IS_ON('WdRHO_P ',myThid) )
142 & doDiagsRho = doDiagsRho + 4
143 IF ( DIAGNOSTICS_IS_ON('WdRHOdP ',myThid) )
144 & doDiagsRho = doDiagsRho + 8
145 ENDIF
146 #endif /* ALLOW_DIAGNOSTICS */
147
148 #ifdef ALLOW_OBCS
149 IF (useOBCS) THEN
150 C-- Calculate future values on open boundaries
151 C-- moved before SEAICE_MODEL call since SEAICE_MODEL needs seaice-obcs fields
152 # ifdef ALLOW_AUTODIFF_TAMC
153 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
154 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
155 # endif
156 # ifdef ALLOW_DEBUG
157 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC',myThid)
158 # endif
159 CALL OBCS_CALC( myTime+deltaTClock, myIter+1,
160 I uVel, vVel, wVel, theta, salt, myThid )
161 ENDIF
162 #endif /* ALLOW_OBCS */
163
164 #ifdef ALLOW_AUTODIFF_TAMC
165 # ifdef ALLOW_SALT_PLUME
166 DO bj=myByLo(myThid),myByHi(myThid)
167 DO bi=myBxLo(myThid),myBxHi(myThid)
168 DO j=1-OLy,sNy+OLy
169 DO i=1-OLx,sNx+OLx
170 saltPlumeDepth(i,j,bi,bj) = 0. _d 0
171 saltPlumeFlux(i,j,bi,bj) = 0. _d 0
172 ENDDO
173 ENDDO
174 ENDDO
175 ENDDO
176 # endif
177 #endif /* ALLOW_AUTODIFF_TAMC */
178
179 #ifdef ALLOW_FRAZIL
180 IF ( useFRAZIL ) THEN
181 C-- Freeze water in the ocean interior and let it rise to the surface
182 CALL FRAZIL_CALC_RHS( myTime, myIter, myThid )
183 ENDIF
184 #endif /* ALLOW_FRAZIL */
185
186 #ifndef OLD_THSICE_CALL_SEQUENCE
187 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
188 IF ( useThSIce .AND. fluidIsWater ) THEN
189 # ifdef ALLOW_AUTODIFF_TAMC
190 CADJ STORE uice,vice = comlev1, key = ikey_dynamics,
191 CADJ & kind = isbyte
192 CADJ STORE iceMask,iceHeight = comlev1, key = ikey_dynamics,
193 CADJ & kind = isbyte
194 CADJ STORE snowHeight, Tsrf = comlev1, key = ikey_dynamics,
195 CADJ & kind = isbyte
196 CADJ STORE Qice1, Qice2 = comlev1, key = ikey_dynamics,
197 CADJ & kind = isbyte
198 CADJ STORE sHeating, snowAge = comlev1, key = ikey_dynamics,
199 CADJ & kind = isbyte
200 CADJ STORE hocemxl = comlev1, key = ikey_dynamics,
201 CADJ & kind = isbyte
202 CADJ STORE icflxsw, snowprc = comlev1, key = ikey_dynamics,
203 CADJ & kind = isbyte
204 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
205 CADJ & kind = isbyte
206 CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics,
207 CADJ & kind = isbyte
208 CADJ STORE qnet,qsw, empmr = comlev1, key = ikey_dynamics,
209 CADJ & kind = isbyte
210 CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics,
211 CADJ & kind = isbyte
212 CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics,
213 CADJ & kind = isbyte
214 # ifdef NONLIN_FRSURF
215 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
216 CADJ & kind = isbyte
217 # endif
218 # endif /* ALLOW_AUTODIFF_TAMC */
219 # ifdef ALLOW_DEBUG
220 IF (debugMode) CALL DEBUG_CALL('THSICE_MAIN',myThid)
221 # endif
222 C-- Step forward Therm.Sea-Ice variables
223 C and modify forcing terms including effects from ice
224 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
225 CALL THSICE_MAIN( myTime, myIter, myThid )
226 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
227 ENDIF
228 #endif /* ALLOW_THSICE */
229 #endif /* ndef OLD_THSICE_CALL_SEQUENCE */
230
231 #ifdef ALLOW_SEAICE
232 # ifdef ALLOW_AUTODIFF_TAMC
233 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
234 CADJ STORE fu,fv = comlev1, key=ikey_dynamics, kind=isbyte
235 CADJ STORE qnet = comlev1, key=ikey_dynamics, kind=isbyte
236 CADJ STORE qsw = comlev1, key=ikey_dynamics, kind=isbyte
237 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
238 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
239 #if (defined ALLOW_EXF) && (defined ALLOW_ATM_TEMP)
240 CADJ STORE evap = comlev1, key=ikey_dynamics, kind=isbyte
241 #endif
242 IF ( .NOT.useSEAICE ) THEN
243 IF ( SEAICEadjMODE .EQ. -1 ) THEN
244 CALL SEAICE_FAKE( myTime, myIter, myThid )
245 ENDIF
246 ENDIF
247 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
248 CADJ STORE fu,fv = comlev1, key=ikey_dynamics, kind=isbyte
249 CADJ STORE qnet = comlev1, key=ikey_dynamics, kind=isbyte
250 CADJ STORE qsw = comlev1, key=ikey_dynamics, kind=isbyte
251 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
252 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
253 #if (defined ALLOW_EXF) && (defined ALLOW_ATM_TEMP)
254 CADJ STORE evap = comlev1, key=ikey_dynamics, kind=isbyte
255 #endif
256 # endif /* ALLOW_AUTODIFF_TAMC */
257 #endif /* ALLOW_SEAICE */
258
259 #ifdef ALLOW_SEAICE
260 IF ( useSEAICE ) THEN
261 # ifdef ALLOW_AUTODIFF_TAMC
262 cph-adj-test(
263 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
264 CADJ STORE hsnow = comlev1, key=ikey_dynamics, kind=isbyte
265 CADJ STORE heff = comlev1, key=ikey_dynamics, kind=isbyte
266 CADJ STORE empmr,qsw,theta = comlev1, key = ikey_dynamics,
267 CADJ & kind = isbyte
268 cph-adj-test)
269 CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics,
270 CADJ & kind = isbyte
271 CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics,
272 CADJ & kind = isbyte
273 cph# ifdef EXF_READ_EVAP
274 CADJ STORE evap = comlev1, key = ikey_dynamics,
275 CADJ & kind = isbyte
276 cph# endif
277 CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics,
278 CADJ & kind = isbyte
279 # ifdef SEAICE_CGRID
280 CADJ STORE stressdivergencex = comlev1, key = ikey_dynamics,
281 CADJ & kind = isbyte
282 CADJ STORE stressdivergencey = comlev1, key = ikey_dynamics,
283 CADJ & kind = isbyte
284 # endif
285 # ifdef SEAICE_ALLOW_DYNAMICS
286 CADJ STORE uice = comlev1, key = ikey_dynamics,
287 CADJ & kind = isbyte
288 CADJ STORE vice = comlev1, key = ikey_dynamics,
289 CADJ & kind = isbyte
290 # ifdef SEAICE_ALLOW_EVP
291 CADJ STORE seaice_sigma1 = comlev1, key = ikey_dynamics,
292 CADJ & kind = isbyte
293 CADJ STORE seaice_sigma2 = comlev1, key = ikey_dynamics,
294 CADJ & kind = isbyte
295 CADJ STORE seaice_sigma12 = comlev1, key = ikey_dynamics,
296 CADJ & kind = isbyte
297 # endif
298 # endif
299 cph# ifdef SEAICE_SALINITY
300 CADJ STORE salt = comlev1, key = ikey_dynamics,
301 CADJ & kind = isbyte
302 cph# endif
303 # ifdef ATMOSPHERIC_LOADING
304 CADJ STORE pload = comlev1, key = ikey_dynamics,
305 CADJ & kind = isbyte
306 CADJ STORE siceload = comlev1, key = ikey_dynamics,
307 CADJ & kind = isbyte
308 # endif
309 # ifdef NONLIN_FRSURF
310 CADJ STORE recip_hfacc = comlev1, key = ikey_dynamics,
311 CADJ & kind = isbyte
312 # endif
313 # ifdef ANNUAL_BALANCE
314 CADJ STORE balance_itcount = comlev1, key = ikey_dynamics,
315 CADJ & kind = isbyte
316 # endif /* ANNUAL_BALANCE */
317 # ifdef ALLOW_THSICE
318 CADJ STORE fu, fv = comlev1, key = ikey_dynamics,
319 CADJ & kind = isbyte
320 # endif
321 # endif /* ALLOW_AUTODIFF_TAMC */
322 # ifdef ALLOW_DEBUG
323 IF (debugMode) CALL DEBUG_CALL('SEAICE_MODEL',myThid)
324 # endif
325 CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
326 CALL SEAICE_MODEL( myTime, myIter, myThid )
327 CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
328 # ifdef ALLOW_COST
329 CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
330 # endif
331 ENDIF
332 #endif /* ALLOW_SEAICE */
333
334 #ifdef ALLOW_AUTODIFF_TAMC
335 CADJ STORE sst, sss = comlev1, key = ikey_dynamics,
336 CADJ & kind = isbyte
337 CADJ STORE qsw = comlev1, key = ikey_dynamics,
338 CADJ & kind = isbyte
339 # ifdef ALLOW_SEAICE
340 CADJ STORE area = comlev1, key = ikey_dynamics,
341 CADJ & kind = isbyte
342 # endif
343 #endif
344
345 #ifdef OLD_THSICE_CALL_SEQUENCE
346 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
347 IF ( useThSIce .AND. fluidIsWater ) THEN
348 # ifdef ALLOW_AUTODIFF_TAMC
349 cph(
350 # ifdef NONLIN_FRSURF
351 CADJ STORE uice,vice = comlev1, key = ikey_dynamics,
352 CADJ & kind = isbyte
353 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
354 CADJ & kind = isbyte
355 CADJ STORE qnet,qsw, empmr = comlev1, key = ikey_dynamics,
356 CADJ & kind = isbyte
357 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
358 CADJ & kind = isbyte
359 # endif
360 # endif
361 # ifdef ALLOW_DEBUG
362 IF (debugMode) CALL DEBUG_CALL('THSICE_MAIN',myThid)
363 # endif
364 C-- Step forward Therm.Sea-Ice variables
365 C and modify forcing terms including effects from ice
366 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
367 CALL THSICE_MAIN( myTime, myIter, myThid )
368 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
369 ENDIF
370 #endif /* ALLOW_THSICE */
371 #endif /* OLD_THSICE_CALL_SEQUENCE */
372
373 #ifdef ALLOW_SHELFICE
374 # ifdef ALLOW_AUTODIFF_TAMC
375 CADJ STORE salt, theta = comlev1, key = ikey_dynamics,
376 CADJ & kind = isbyte
377 # endif
378 IF ( useShelfIce .AND. fluidIsWater ) THEN
379 #ifdef ALLOW_DEBUG
380 IF (debugMode) CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
381 #endif
382 C compute temperature and (virtual) salt flux at the
383 C shelf-ice ocean interface
384 CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
385 & myThid)
386 CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
387 CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
388 & myThid)
389 ENDIF
390 #endif /* ALLOW_SHELFICE */
391
392 #ifdef ALLOW_ICEFRONT
393 IF ( useICEFRONT .AND. fluidIsWater ) THEN
394 #ifdef ALLOW_DEBUG
395 IF (debugMode) CALL DEBUG_CALL('ICEFRONT_THERMODYNAMICS',myThid)
396 #endif
397 C compute temperature and (virtual) salt flux at the
398 C ice-front ocean interface
399 CALL TIMER_START('ICEFRONT_THERMODYNAMICS [DO_OCEANIC_PHYS]',
400 & myThid)
401 CALL ICEFRONT_THERMODYNAMICS( myTime, myIter, myThid )
402 CALL TIMER_STOP( 'ICEFRONT_THERMODYNAMICS [DO_OCEANIC_PHYS]',
403 & myThid)
404 ENDIF
405 #endif /* ALLOW_ICEFRONT */
406
407 #ifdef ALLOW_SALT_PLUME
408 IF ( useSALT_PLUME ) THEN
409 CALL SALT_PLUME_DO_EXCH( myTime, myIter, myThid )
410 ENDIF
411 #endif /* ALLOW_SALT_PLUME */
412
413 C-- Freeze water at the surface
414 IF ( allowFreezing ) THEN
415 #ifdef ALLOW_AUTODIFF_TAMC
416 CADJ STORE theta = comlev1, key = ikey_dynamics,
417 CADJ & kind = isbyte
418 #endif
419 CALL FREEZE_SURFACE( myTime, myIter, myThid )
420 ENDIF
421
422 #ifdef ALLOW_OCN_COMPON_INTERF
423 C-- Apply imported data (from coupled interface) to forcing fields
424 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
425 IF ( useCoupler ) THEN
426 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
427 ENDIF
428 #endif /* ALLOW_OCN_COMPON_INTERF */
429
430 #ifdef ALLOW_BALANCE_FLUXES
431 C balance fluxes
432 IF ( balanceEmPmR .AND. (.NOT.useSeaice .OR. useThSIce) )
433 & CALL REMOVE_MEAN_RS( 1, EmPmR, maskInC, maskInC, rA, drF,
434 & 'EmPmR', myTime, myThid )
435 IF ( balanceQnet .AND. (.NOT.useSeaice .OR. useThSIce) )
436 & CALL REMOVE_MEAN_RS( 1, Qnet, maskInC, maskInC, rA, drF,
437 & 'Qnet ', myTime, myThid )
438 #endif /* ALLOW_BALANCE_FLUXES */
439
440 #ifdef ALLOW_AUTODIFF_TAMC
441 C-- HPF directive to help TAMC
442 CHPF$ INDEPENDENT
443 #else /* ALLOW_AUTODIFF_TAMC */
444 C if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
445 C and all vertical mixing schemes, but keep OBCS_CALC
446 IF ( fluidIsWater ) THEN
447 #endif /* ALLOW_AUTODIFF_TAMC */
448 DO bj=myByLo(myThid),myByHi(myThid)
449 #ifdef ALLOW_AUTODIFF_TAMC
450 C-- HPF directive to help TAMC
451 CHPF$ INDEPENDENT
452 #endif /* ALLOW_AUTODIFF_TAMC */
453 DO bi=myBxLo(myThid),myBxHi(myThid)
454
455 #ifdef ALLOW_AUTODIFF_TAMC
456 act1 = bi - myBxLo(myThid)
457 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
458 act2 = bj - myByLo(myThid)
459 max2 = myByHi(myThid) - myByLo(myThid) + 1
460 act3 = myThid - 1
461 max3 = nTx*nTy
462 act4 = ikey_dynamics - 1
463 itdkey = (act1 + 1) + act2*max1
464 & + act3*max1*max2
465 & + act4*max1*max2*max3
466 #endif /* ALLOW_AUTODIFF_TAMC */
467
468 C-- Set up work arrays with valid (i.e. not NaN) values
469 C These inital values do not alter the numerical results. They
470 C just ensure that all memory references are to valid floating
471 C point numbers. This prevents spurious hardware signals due to
472 C uninitialised but inert locations.
473
474 #ifdef ALLOW_AUTODIFF_TAMC
475 DO j=1-OLy,sNy+OLy
476 DO i=1-OLx,sNx+OLx
477 rhoKm1 (i,j) = 0. _d 0
478 rhoKp1 (i,j) = 0. _d 0
479 ENDDO
480 ENDDO
481 #endif /* ALLOW_AUTODIFF_TAMC */
482
483 DO k=1,Nr
484 DO j=1-OLy,sNy+OLy
485 DO i=1-OLx,sNx+OLx
486 C This is currently used by GMRedi, IVDC, MXL-depth and Diagnostics
487 sigmaX(i,j,k) = 0. _d 0
488 sigmaY(i,j,k) = 0. _d 0
489 sigmaR(i,j,k) = 0. _d 0
490 #ifdef ALLOW_AUTODIFF_TAMC
491 cph all the following init. are necessary for TAF
492 cph although some of these are re-initialised later.
493 rhoInSitu(i,j,k,bi,bj) = 0.
494 IVDConvCount(i,j,k,bi,bj) = 0.
495 # ifdef ALLOW_GMREDI
496 Kwx(i,j,k,bi,bj) = 0. _d 0
497 Kwy(i,j,k,bi,bj) = 0. _d 0
498 Kwz(i,j,k,bi,bj) = 0. _d 0
499 # ifdef GM_NON_UNITY_DIAGONAL
500 Kux(i,j,k,bi,bj) = 0. _d 0
501 Kvy(i,j,k,bi,bj) = 0. _d 0
502 # endif
503 # ifdef GM_EXTRA_DIAGONAL
504 Kuz(i,j,k,bi,bj) = 0. _d 0
505 Kvz(i,j,k,bi,bj) = 0. _d 0
506 # endif
507 # ifdef GM_BOLUS_ADVEC
508 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
509 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
510 # endif
511 # ifdef GM_VISBECK_VARIABLE_K
512 VisbeckK(i,j,bi,bj) = 0. _d 0
513 # endif
514 # endif /* ALLOW_GMREDI */
515 # ifdef ALLOW_KPP
516 KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0
517 KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0
518 # endif /* ALLOW_KPP */
519 # ifdef ALLOW_GGL90
520 GGL90viscArU(i,j,k,bi,bj) = 0. _d 0
521 GGL90viscArV(i,j,k,bi,bj) = 0. _d 0
522 GGL90diffKr(i,j,k,bi,bj) = 0. _d 0
523 # endif /* ALLOW_GGL90 */
524 #endif /* ALLOW_AUTODIFF_TAMC */
525 ENDDO
526 ENDDO
527 ENDDO
528
529 iMin = 1-OLx
530 iMax = sNx+OLx
531 jMin = 1-OLy
532 jMax = sNy+OLy
533
534 #ifdef ALLOW_AUTODIFF_TAMC
535 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
536 CADJ & kind = isbyte
537 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
538 CADJ & kind = isbyte
539 CADJ STORE totphihyd(:,:,:,bi,bj)
540 CADJ & = comlev1_bibj, key=itdkey,
541 CADJ & kind = isbyte
542 # ifdef ALLOW_KPP
543 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
544 CADJ & kind = isbyte
545 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
546 CADJ & kind = isbyte
547 # endif
548 # ifdef ALLOW_SALT_PLUME
549 CADJ STORE saltplumedepth(:,:,bi,bj) = comlev1_bibj, key=itdkey,
550 CADJ & kind = isbyte
551 # endif
552 #endif /* ALLOW_AUTODIFF_TAMC */
553
554 C-- Always compute density (stored in common block) here; even when it is not
555 C needed here, will be used anyway in calc_phi_hyd (data flow easier this way)
556 #ifdef ALLOW_DEBUG
557 IF (debugMode) CALL DEBUG_CALL('FIND_RHO_2D (xNr)',myThid)
558 #endif
559 #ifdef ALLOW_AUTODIFF_TAMC
560 IF ( fluidIsWater ) THEN
561 #endif /* ALLOW_AUTODIFF_TAMC */
562 #ifdef ALLOW_DOWN_SLOPE
563 IF ( useDOWN_SLOPE ) THEN
564 DO k=1,Nr
565 CALL DWNSLP_CALC_RHO(
566 I theta, salt,
567 O rhoInSitu(1-OLx,1-OLy,k,bi,bj),
568 I k, bi, bj, myTime, myIter, myThid )
569 ENDDO
570 ENDIF
571 #endif /* ALLOW_DOWN_SLOPE */
572 #ifdef ALLOW_BBL
573 IF ( useBBL ) THEN
574 C pkg/bbl requires in-situ bbl density for depths equal to and deeper than the bbl.
575 C To reduce computation and storage requirement, these densities are stored in the
576 C dry grid boxes of rhoInSitu. See BBL_CALC_RHO for details.
577 DO k=Nr,1,-1
578 CALL BBL_CALC_RHO(
579 I theta, salt,
580 O rhoInSitu,
581 I k, bi, bj, myTime, myIter, myThid )
582
583 ENDDO
584 ENDIF
585 #endif /* ALLOW_BBL */
586 IF ( .NOT. ( useDOWN_SLOPE .OR. useBBL ) ) THEN
587 DO k=1,Nr
588 CALL FIND_RHO_2D(
589 I iMin, iMax, jMin, jMax, k,
590 I theta(1-OLx,1-OLy,k,bi,bj),
591 I salt (1-OLx,1-OLy,k,bi,bj),
592 O rhoInSitu(1-OLx,1-OLy,k,bi,bj),
593 I k, bi, bj, myThid )
594 ENDDO
595 ENDIF
596 #ifdef ALLOW_AUTODIFF_TAMC
597 ELSE
598 C- fluid is not water:
599 DO k=1,Nr
600 DO j=1-OLy,sNy+OLy
601 DO i=1-OLx,sNx+OLx
602 rhoInSitu(i,j,k,bi,bj) = 0.
603 ENDDO
604 ENDDO
605 ENDDO
606 ENDIF
607 #endif /* ALLOW_AUTODIFF_TAMC */
608
609 #ifdef ALLOW_DEBUG
610 IF (debugMode) CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
611 #endif
612
613 C-- Start of diagnostic loop
614 DO k=Nr,1,-1
615
616 #ifdef ALLOW_AUTODIFF_TAMC
617 C? Patrick, is this formula correct now that we change the loop range?
618 C? Do we still need this?
619 cph kkey formula corrected.
620 cph Needed for rhoK, rhoKm1, in the case useGMREDI.
621 kkey = (itdkey-1)*Nr + k
622 #endif /* ALLOW_AUTODIFF_TAMC */
623
624 c#ifdef ALLOW_AUTODIFF_TAMC
625 cCADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
626 cCADJ & kind = isbyte
627 cCADJ STORE salt(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
628 cCADJ & kind = isbyte
629 c#endif /* ALLOW_AUTODIFF_TAMC */
630
631 C-- Calculate gradients of potential density for isoneutral
632 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
633 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
634 & .OR. usePP81 .OR. useMY82 .OR. useGGL90
635 & .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
636 IF (k.GT.1) THEN
637 #ifdef ALLOW_AUTODIFF_TAMC
638 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey,
639 CADJ & kind = isbyte
640 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey,
641 CADJ & kind = isbyte
642 CADJ STORE rhokm1 (bi,bj) = comlev1_bibj_k, key=kkey,
643 CADJ & kind = isbyte
644 #endif /* ALLOW_AUTODIFF_TAMC */
645 CALL FIND_RHO_2D(
646 I iMin, iMax, jMin, jMax, k,
647 I theta(1-OLx,1-OLy,k-1,bi,bj),
648 I salt (1-OLx,1-OLy,k-1,bi,bj),
649 O rhoKm1,
650 I k-1, bi, bj, myThid )
651 ENDIF
652 #ifdef ALLOW_DEBUG
653 IF (debugMode) CALL DEBUG_CALL('GRAD_SIGMA',myThid)
654 #endif
655 cph Avoid variable aliasing for adjoint !!!
656 DO j=jMin,jMax
657 DO i=iMin,iMax
658 rhoKp1(i,j) = rhoInSitu(i,j,k,bi,bj)
659 ENDDO
660 ENDDO
661 CALL GRAD_SIGMA(
662 I bi, bj, iMin, iMax, jMin, jMax, k,
663 I rhoInSitu(1-OLx,1-OLy,k,bi,bj), rhoKm1, rhoKp1,
664 O sigmaX, sigmaY, sigmaR,
665 I myThid )
666 #ifdef ALLOW_AUTODIFF_TAMC
667 #ifdef GMREDI_WITH_STABLE_ADJOINT
668 cgf zero out adjoint fields to stabilize pkg/gmredi adjoint
669 cgf -> cuts adjoint dependency from slope to state
670 CALL ZERO_ADJ_LOC( Nr, sigmaX, myThid)
671 CALL ZERO_ADJ_LOC( Nr, sigmaY, myThid)
672 CALL ZERO_ADJ_LOC( Nr, sigmaR, myThid)
673 #endif
674 #endif /* ALLOW_AUTODIFF_TAMC */
675 ENDIF
676
677 C-- Implicit Vertical Diffusion for Convection
678 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
679 #ifdef ALLOW_DEBUG
680 IF (debugMode) CALL DEBUG_CALL('CALC_IVDC',myThid)
681 #endif
682 CALL CALC_IVDC(
683 I bi, bj, iMin, iMax, jMin, jMax, k,
684 I sigmaR,
685 I myTime, myIter, myThid)
686 ENDIF
687
688 #ifdef ALLOW_DIAGNOSTICS
689 IF ( doDiagsRho.GE.4 ) THEN
690 CALL DIAGS_RHO_L( doDiagsRho, k, bi, bj,
691 I rhoInSitu(1-OLx,1-OLy,1,bi,bj),
692 I rhoKm1, wVel,
693 I myTime, myIter, myThid )
694 ENDIF
695 #endif
696
697 C-- end of diagnostic k loop (Nr:1)
698 ENDDO
699
700 #ifdef ALLOW_AUTODIFF_TAMC
701 CADJ STORE IVDConvCount(:,:,:,bi,bj)
702 CADJ & = comlev1_bibj, key=itdkey,
703 CADJ & kind = isbyte
704 #endif
705
706 C-- Diagnose Mixed Layer Depth:
707 IF ( useGMRedi .OR. MOD(doDiagsRho,2).EQ.1 ) THEN
708 CALL CALC_OCE_MXLAYER(
709 I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR,
710 I bi, bj, myTime, myIter, myThid )
711 ENDIF
712
713 #ifdef ALLOW_SALT_PLUME
714 IF ( useSALT_PLUME ) THEN
715 CALL SALT_PLUME_CALC_DEPTH(
716 I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR,
717 I bi, bj, myTime, myIter, myThid )
718 ENDIF
719 #endif /* ALLOW_SALT_PLUME */
720
721 #ifdef ALLOW_DIAGNOSTICS
722 IF ( MOD(doDiagsRho,4).GE.2 ) THEN
723 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
724 & 2, bi, bj, myThid)
725 ENDIF
726 #endif /* ALLOW_DIAGNOSTICS */
727
728 C-- Determines forcing terms based on external fields
729 C relaxation terms, etc.
730 #ifdef ALLOW_DEBUG
731 IF (debugMode) CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
732 #endif
733 #ifdef ALLOW_AUTODIFF_TAMC
734 CADJ STORE EmPmR(:,:,bi,bj)
735 CADJ & = comlev1_bibj, key=itdkey,
736 CADJ & kind = isbyte
737 # ifdef EXACT_CONSERV
738 CADJ STORE PmEpR(:,:,bi,bj)
739 CADJ & = comlev1_bibj, key=itdkey,
740 CADJ & kind = isbyte
741 # endif
742 # ifdef NONLIN_FRSURF
743 CADJ STORE hFac_surfC(:,:,bi,bj)
744 CADJ & = comlev1_bibj, key=itdkey,
745 CADJ & kind = isbyte
746 CADJ STORE recip_hFacC(:,:,:,bi,bj)
747 CADJ & = comlev1_bibj, key=itdkey,
748 CADJ & kind = isbyte
749 # if (defined (ALLOW_PTRACERS))
750 CADJ STORE surfaceForcingS(:,:,bi,bj) = comlev1_bibj, key=itdkey,
751 CADJ & kind = isbyte
752 CADJ STORE surfaceForcingT(:,:,bi,bj) = comlev1_bibj, key=itdkey,
753 CADJ & kind = isbyte
754 # endif /* ALLOW_PTRACERS */
755 # endif /* NONLIN_FRSURF */
756 #endif
757 CALL EXTERNAL_FORCING_SURF(
758 I bi, bj, iMin, iMax, jMin, jMax,
759 I myTime, myIter, myThid )
760 #ifdef ALLOW_AUTODIFF_TAMC
761 # ifdef EXACT_CONSERV
762 cph-test
763 cphCADJ STORE PmEpR(:,:,bi,bj)
764 cphCADJ & = comlev1_bibj, key=itdkey,
765 cphCADJ & kind = isbyte
766 # endif
767 #endif
768
769 #ifdef ALLOW_AUTODIFF_TAMC
770 cph needed for KPP
771 CADJ STORE surfaceForcingU(:,:,bi,bj)
772 CADJ & = comlev1_bibj, key=itdkey,
773 CADJ & kind = isbyte
774 CADJ STORE surfaceForcingV(:,:,bi,bj)
775 CADJ & = comlev1_bibj, key=itdkey,
776 CADJ & kind = isbyte
777 CADJ STORE surfaceForcingS(:,:,bi,bj)
778 CADJ & = comlev1_bibj, key=itdkey,
779 CADJ & kind = isbyte
780 CADJ STORE surfaceForcingT(:,:,bi,bj)
781 CADJ & = comlev1_bibj, key=itdkey,
782 CADJ & kind = isbyte
783 CADJ STORE surfaceForcingTice(:,:,bi,bj)
784 CADJ & = comlev1_bibj, key=itdkey,
785 CADJ & kind = isbyte
786 #endif /* ALLOW_AUTODIFF_TAMC */
787
788 #ifdef ALLOW_KPP
789 C-- Compute KPP mixing coefficients
790 IF (useKPP) THEN
791 #ifdef ALLOW_DEBUG
792 IF (debugMode) CALL DEBUG_CALL('KPP_CALC',myThid)
793 #endif
794 CALL TIMER_START('KPP_CALC [DO_OCEANIC_PHYS]', myThid)
795 CALL KPP_CALC(
796 I bi, bj, myTime, myIter, myThid )
797 CALL TIMER_STOP ('KPP_CALC [DO_OCEANIC_PHYS]', myThid)
798 #ifdef ALLOW_AUTODIFF_TAMC
799 ELSE
800 CALL KPP_CALC_DUMMY(
801 I bi, bj, myTime, myIter, myThid )
802 #endif /* ALLOW_AUTODIFF_TAMC */
803 ENDIF
804
805 #endif /* ALLOW_KPP */
806
807 #ifdef ALLOW_PP81
808 C-- Compute PP81 mixing coefficients
809 IF (usePP81) THEN
810 #ifdef ALLOW_DEBUG
811 IF (debugMode) CALL DEBUG_CALL('PP81_CALC',myThid)
812 #endif
813 CALL PP81_CALC(
814 I bi, bj, sigmaR, myTime, myIter, myThid )
815 ENDIF
816 #endif /* ALLOW_PP81 */
817
818 #ifdef ALLOW_MY82
819 C-- Compute MY82 mixing coefficients
820 IF (useMY82) THEN
821 #ifdef ALLOW_DEBUG
822 IF (debugMode) CALL DEBUG_CALL('MY82_CALC',myThid)
823 #endif
824 CALL MY82_CALC(
825 I bi, bj, sigmaR, myTime, myIter, myThid )
826 ENDIF
827 #endif /* ALLOW_MY82 */
828
829 #ifdef ALLOW_GGL90
830 #ifdef ALLOW_AUTODIFF_TAMC
831 CADJ STORE GGL90TKE (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
832 CADJ & kind = isbyte
833 #endif /* ALLOW_AUTODIFF_TAMC */
834 C-- Compute GGL90 mixing coefficients
835 IF (useGGL90) THEN
836 #ifdef ALLOW_DEBUG
837 IF (debugMode) CALL DEBUG_CALL('GGL90_CALC',myThid)
838 #endif
839 CALL TIMER_START('GGL90_CALC [DO_OCEANIC_PHYS]', myThid)
840 CALL GGL90_CALC(
841 I bi, bj, sigmaR, myTime, myIter, myThid )
842 CALL TIMER_STOP ('GGL90_CALC [DO_OCEANIC_PHYS]', myThid)
843 ENDIF
844 #endif /* ALLOW_GGL90 */
845
846 #ifdef ALLOW_TIMEAVE
847 IF ( taveFreq.GT. 0. _d 0 ) THEN
848 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
849 ENDIF
850 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
851 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
852 I Nr, deltaTClock, bi, bj, myThid)
853 ENDIF
854 #endif /* ALLOW_TIMEAVE */
855
856 #ifdef ALLOW_GMREDI
857 #ifdef ALLOW_AUTODIFF_TAMC
858 # ifndef GM_EXCLUDE_CLIPPING
859 cph storing here is needed only for one GMREDI_OPTIONS:
860 cph define GM_BOLUS_ADVEC
861 cph keep it although TAF says you dont need to.
862 cph but I have avoided the #ifdef for now, in case more things change
863 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey,
864 CADJ & kind = isbyte
865 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey,
866 CADJ & kind = isbyte
867 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey,
868 CADJ & kind = isbyte
869 # endif
870 #endif /* ALLOW_AUTODIFF_TAMC */
871
872 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
873 IF (useGMRedi) THEN
874 #ifdef ALLOW_DEBUG
875 IF (debugMode) CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
876 #endif
877 CALL GMREDI_CALC_TENSOR(
878 I iMin, iMax, jMin, jMax,
879 I sigmaX, sigmaY, sigmaR,
880 I bi, bj, myTime, myIter, myThid )
881 #ifdef ALLOW_AUTODIFF_TAMC
882 ELSE
883 CALL GMREDI_CALC_TENSOR_DUMMY(
884 I iMin, iMax, jMin, jMax,
885 I sigmaX, sigmaY, sigmaR,
886 I bi, bj, myTime, myIter, myThid )
887 #endif /* ALLOW_AUTODIFF_TAMC */
888 ENDIF
889 #endif /* ALLOW_GMREDI */
890
891 #ifdef ALLOW_DOWN_SLOPE
892 IF ( useDOWN_SLOPE ) THEN
893 C-- Calculate Downsloping Flow for Down_Slope parameterization
894 IF ( usingPCoords ) THEN
895 CALL DWNSLP_CALC_FLOW(
896 I bi, bj, kSurfC, rhoInSitu,
897 I myTime, myIter, myThid )
898 ELSE
899 CALL DWNSLP_CALC_FLOW(
900 I bi, bj, kLowC, rhoInSitu,
901 I myTime, myIter, myThid )
902 ENDIF
903 ENDIF
904 #endif /* ALLOW_DOWN_SLOPE */
905
906 C-- end bi,bj loops.
907 ENDDO
908 ENDDO
909
910 #ifdef ALLOW_BALANCE_RELAX
911 # ifdef ALLOW_AUTODIFF_TAMC
912 CADJ STORE SSSrlx = comlev1, key=ikey_dynamics, kind=isbyte
913 CADJ STORE SSSrlxTile = comlev1, key=ikey_dynamics, kind=isbyte
914 CADJ STORE SSSrlxGlob = comlev1, key=ikey_dynamics, kind=isbyte
915 CADJ STORE SSTrlx = comlev1, key=ikey_dynamics, kind=isbyte
916 CADJ STORE SSTrlxTile = comlev1, key=ikey_dynamics, kind=isbyte
917 CADJ STORE SSTrlxGlob = comlev1, key=ikey_dynamics, kind=isbyte
918 # endif /* ALLOW_AUTODIFF_TAMC */
919 CALL BALANCE_RELAX( myTime, myIter, myThid )
920 #endif /* ALLOW_BALANCE_RELAX */
921
922 #ifndef ALLOW_AUTODIFF_TAMC
923 C--- if fluid Is Water: end
924 ENDIF
925 #endif
926
927 #ifdef ALLOW_BBL
928 IF ( useBBL ) THEN
929 CALL BBL_CALC_RHS(
930 I myTime, myIter, myThid )
931 ENDIF
932 #endif /* ALLOW_BBL */
933
934 #ifdef ALLOW_MYPACKAGE
935 IF ( useMYPACKAGE ) THEN
936 CALL MYPACKAGE_CALC_RHS(
937 I myTime, myIter, myThid )
938 ENDIF
939 #endif /* ALLOW_MYPACKAGE */
940
941 #ifdef ALLOW_GMREDI
942 IF ( useGMRedi ) THEN
943 CALL GMREDI_DO_EXCH( myTime, myIter, myThid )
944 ENDIF
945 #endif /* ALLOW_GMREDI */
946
947 #ifdef ALLOW_KPP
948 IF (useKPP) THEN
949 CALL KPP_DO_EXCH( myThid )
950 ENDIF
951 #endif /* ALLOW_KPP */
952
953 #ifdef ALLOW_DIAGNOSTICS
954 IF ( fluidIsWater .AND. useDiagnostics ) THEN
955 CALL DIAGS_RHO_G(
956 I rhoInSitu, uVel, vVel, wVel,
957 I myTime, myIter, myThid )
958 CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
959 ENDIF
960 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
961 CALL DIAGNOSTICS_FILL( IVDConvCount, 'CONVADJ ',
962 & 0, Nr, 0, 1, 1, myThid )
963 ENDIF
964 #endif
965
966 #ifdef ALLOW_ECCO
967 CALL ECCO_PHYS(mythid)
968 #endif
969
970 #ifdef ALLOW_DEBUG
971 IF (debugMode) CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
972 #endif
973
974 RETURN
975 END

  ViewVC Help
Powered by ViewVC 1.1.22