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

Annotation of /MITgcm/model/src/dynamics.F

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


Revision 1.49 - (hide annotations) (download)
Fri Jun 9 02:45:04 2000 UTC (23 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint28
Changes since 1.48: +253 -11 lines
Modifications to include TAMC directives, tape key computations
and initialisations to make code TAMC compatible.
Routines the_model_main.F and initialise_fixed.F
are left unchanged for the moment. (P.H.)

1 heimbach 1.49 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.48 2000/03/14 17:47:25 adcroft Exp $
2 cnh 1.1
3 adcroft 1.24 #include "CPP_OPTIONS.h"
4 cnh 1.1
5 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6 cnh 1.1 C /==========================================================\
7     C | SUBROUTINE DYNAMICS |
8     C | o Controlling routine for the explicit part of the model |
9     C | dynamics. |
10     C |==========================================================|
11     C | This routine evaluates the "dynamics" terms for each |
12     C | block of ocean in turn. Because the blocks of ocean have |
13     C | overlap regions they are independent of one another. |
14     C | If terms involving lateral integrals are needed in this |
15     C | routine care will be needed. Similarly finite-difference |
16     C | operations with stencils wider than the overlap region |
17     C | require special consideration. |
18     C | Notes |
19     C | ===== |
20     C | C*P* comments indicating place holders for which code is |
21     C | presently being developed. |
22     C \==========================================================/
23 heimbach 1.49 c
24     c changed: Patrick Heimbach heimbach@mit.edu 6-Jun-2000
25     c - computation of ikey wrong for nTx,nTy > 1
26     c and/or nsx,nsy > 1: act1 and act2 were
27     c mixed up.
28    
29 adcroft 1.40 IMPLICIT NONE
30 cnh 1.1
31     C == Global variables ===
32     #include "SIZE.h"
33     #include "EEPARAMS.h"
34     #include "CG2D.h"
35 adcroft 1.6 #include "PARAMS.h"
36 adcroft 1.3 #include "DYNVARS.h"
37 adcroft 1.42 #include "GRID.h"
38 heimbach 1.49
39 adcroft 1.42 #ifdef ALLOW_KPP
40     #include "KPPMIX.h"
41     #endif
42 cnh 1.1
43 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
44     #include "tamc.h"
45     #include "tamc_keys.h"
46     #endif
47    
48 cnh 1.1 C == Routine arguments ==
49 cnh 1.8 C myTime - Current time in simulation
50     C myIter - Current iteration number in simulation
51 cnh 1.1 C myThid - Thread number for this instance of the routine.
52 cnh 1.8 _RL myTime
53     INTEGER myIter
54 adcroft 1.47 INTEGER myThid
55 cnh 1.1
56     C == Local variables
57     C xA, yA - Per block temporaries holding face areas
58 cnh 1.38 C uTrans, vTrans, rTrans - Per block temporaries holding flow
59     C transport
60 cnh 1.30 C rVel o uTrans: Zonal transport
61 cnh 1.1 C o vTrans: Meridional transport
62 cnh 1.30 C o rTrans: Vertical transport
63 cnh 1.38 C o rVel: Vertical velocity at upper and
64     C lower cell faces.
65 cnh 1.1 C maskC,maskUp o maskC: land/water mask for tracer cells
66     C o maskUp: land/water mask for W points
67     C aTerm, xTerm, cTerm - Work arrays for holding separate terms in
68     C mTerm, pTerm, tendency equations.
69     C fZon, fMer, fVer[STUV] o aTerm: Advection term
70     C o xTerm: Mixing term
71     C o cTerm: Coriolis term
72     C o mTerm: Metric term
73     C o pTerm: Pressure term
74     C o fZon: Zonal flux term
75     C o fMer: Meridional flux term
76     C o fVer: Vertical flux term - note fVer
77     C is "pipelined" in the vertical
78     C so we need an fVer for each
79     C variable.
80 cnh 1.38 C rhoK, rhoKM1 - Density at current level, level above and level
81     C below.
82 cnh 1.26 C rhoKP1
83     C buoyK, buoyKM1 - Buoyancy at current level and level above.
84 cnh 1.31 C phiHyd - Hydrostatic part of the potential phiHydi.
85 cnh 1.38 C In z coords phiHydiHyd is the hydrostatic
86     C pressure anomaly
87     C In p coords phiHydiHyd is the geopotential
88     C surface height
89 cnh 1.30 C anomaly.
90     C etaSurfX, - Holds surface elevation gradient in X and Y.
91     C etaSurfY
92     C K13, K23, K33 - Non-zero elements of small-angle approximation
93     C diffusion tensor.
94     C KapGM - Spatially varying Visbeck et. al mixing coeff.
95     C KappaRT, - Total diffusion in vertical for T and S.
96 cnh 1.38 C KappaRS (background + spatially varying, isopycnal term).
97 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
98     C jMin, jMax are applied.
99 cnh 1.1 C bi, bj
100 cnh 1.30 C k, kUp, - Index for layer above and below. kUp and kDown
101 cnh 1.38 C kDown, kM1 are switched with layer to be the appropriate
102     C index into fVerTerm.
103 cnh 1.30 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108     _RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
109     _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
111     _RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112     _RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113     _RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114     _RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115     _RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118     _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
119     _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
120     _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
121     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
122 cnh 1.31 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
123 cnh 1.30 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124     _RL rhokp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126     _RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127     _RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128     _RL rhotmp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 cnh 1.29 _RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130     _RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 cnh 1.31 _RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
132     _RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
133     _RL K33 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
134 cnh 1.30 _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135 cnh 1.31 _RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
136     _RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
137 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
138     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
139 adcroft 1.12
140 adcroft 1.45 #ifdef INCLUDE_CONVECT_CALL
141     _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
142     #endif
143    
144 cnh 1.1 INTEGER iMin, iMax
145     INTEGER jMin, jMax
146     INTEGER bi, bj
147     INTEGER i, j
148     INTEGER k, kM1, kUp, kDown
149 cnh 1.19 LOGICAL BOTTOM_LAYER
150 cnh 1.1
151 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
152     INTEGER isbyte
153     PARAMETER( isbyte = 4 )
154    
155     INTEGER act1, act2, act3, act4
156     INTEGER max1, max2, max3
157     INTEGER ikact, iikey,kkey
158     INTEGER maximpl
159     #endif
160    
161 adcroft 1.11 C--- The algorithm...
162     C
163     C "Correction Step"
164     C =================
165     C Here we update the horizontal velocities with the surface
166     C pressure such that the resulting flow is either consistent
167     C with the free-surface evolution or the rigid-lid:
168     C U[n] = U* + dt x d/dx P
169     C V[n] = V* + dt x d/dy P
170     C
171     C "Calculation of Gs"
172     C ===================
173     C This is where all the accelerations and tendencies (ie.
174 cnh 1.31 C phiHydysics, parameterizations etc...) are calculated
175 cnh 1.27 C rVel = sum_r ( div. u[n] )
176 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
177 cnh 1.27 C b = b(rho, theta)
178 adcroft 1.11 C K31 = K31 ( rho )
179 cnh 1.27 C Gu[n] = Gu( u[n], v[n], rVel, b, ... )
180     C Gv[n] = Gv( u[n], v[n], rVel, b, ... )
181     C Gt[n] = Gt( theta[n], u[n], v[n], rVel, K31, ... )
182     C Gs[n] = Gs( salt[n], u[n], v[n], rVel, K31, ... )
183 adcroft 1.11 C
184 adcroft 1.12 C "Time-stepping" or "Prediction"
185 adcroft 1.11 C ================================
186     C The models variables are stepped forward with the appropriate
187     C time-stepping scheme (currently we use Adams-Bashforth II)
188     C - For momentum, the result is always *only* a "prediction"
189     C in that the flow may be divergent and will be "corrected"
190     C later with a surface pressure gradient.
191     C - Normally for tracers the result is the new field at time
192     C level [n+1} *BUT* in the case of implicit diffusion the result
193     C is also *only* a prediction.
194     C - We denote "predictors" with an asterisk (*).
195     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
196     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
197     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
198     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
199 adcroft 1.12 C With implicit diffusion:
200 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
201     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
202 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
203     C (1 + dt * K * d_zz) salt[n] = salt*
204 adcroft 1.11 C---
205    
206 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
207     C-- dummy statement to end declaration part
208     ikey = 1
209     #endif
210    
211 cnh 1.1 C-- Set up work arrays with valid (i.e. not NaN) values
212     C These inital values do not alter the numerical results. They
213     C just ensure that all memory references are to valid floating
214     C point numbers. This prevents spurious hardware signals due to
215     C uninitialised but inert locations.
216     DO j=1-OLy,sNy+OLy
217     DO i=1-OLx,sNx+OLx
218 adcroft 1.5 xA(i,j) = 0. _d 0
219     yA(i,j) = 0. _d 0
220     uTrans(i,j) = 0. _d 0
221     vTrans(i,j) = 0. _d 0
222     aTerm(i,j) = 0. _d 0
223     xTerm(i,j) = 0. _d 0
224     cTerm(i,j) = 0. _d 0
225     mTerm(i,j) = 0. _d 0
226     pTerm(i,j) = 0. _d 0
227     fZon(i,j) = 0. _d 0
228     fMer(i,j) = 0. _d 0
229 cnh 1.31 DO K=1,Nr
230     phiHyd (i,j,k) = 0. _d 0
231 cnh 1.30 K13(i,j,k) = 0. _d 0
232     K23(i,j,k) = 0. _d 0
233     K33(i,j,k) = 0. _d 0
234 adcroft 1.45 KappaRU(i,j,k) = 0. _d 0
235     KappaRV(i,j,k) = 0. _d 0
236 cnh 1.1 ENDDO
237 cnh 1.30 rhoKM1 (i,j) = 0. _d 0
238     rhok (i,j) = 0. _d 0
239     rhoKP1 (i,j) = 0. _d 0
240     rhoTMP (i,j) = 0. _d 0
241 cnh 1.26 buoyKM1(i,j) = 0. _d 0
242     buoyK (i,j) = 0. _d 0
243 cnh 1.30 maskC (i,j) = 0. _d 0
244 cnh 1.1 ENDDO
245     ENDDO
246    
247 cnh 1.35
248 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
249     C-- HPF directive to help TAMC
250     !HPF$ INDEPENDENT
251     #endif
252    
253 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
254 heimbach 1.49
255     #ifdef ALLOW_AUTODIFF_TAMC
256     C-- HPF directive to help TAMC
257     !HPF$ INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV
258     !HPF$& ,phiHyd,K13,K23,K33,KapGM
259     !HPF$& ,utrans,vtrans,maskc,xA,yA
260     !HPF$& ,KappaRT,KappaRS,KappaRU,KappaRV
261     !HPF$& )
262     #endif
263    
264 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
265    
266 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
267     act1 = bi - myBxLo(myThid)
268     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
269    
270     act2 = bj - myByLo(myThid)
271     max2 = myByHi(myThid) - myByLo(myThid) + 1
272    
273     act3 = myThid - 1
274     max3 = nTx*nTy
275    
276     act4 = ikey_dynamics - 1
277    
278     ikey = (act1 + 1) + act2*max1
279     & + act3*max1*max2
280     & + act4*max1*max2*max3
281     #endif
282    
283 cnh 1.7 C-- Set up work arrays that need valid initial values
284     DO j=1-OLy,sNy+OLy
285     DO i=1-OLx,sNx+OLx
286 cnh 1.27 rTrans(i,j) = 0. _d 0
287     rVel (i,j,1) = 0. _d 0
288     rVel (i,j,2) = 0. _d 0
289 cnh 1.30 fVerT (i,j,1) = 0. _d 0
290     fVerT (i,j,2) = 0. _d 0
291     fVerS (i,j,1) = 0. _d 0
292     fVerS (i,j,2) = 0. _d 0
293     fVerU (i,j,1) = 0. _d 0
294     fVerU (i,j,2) = 0. _d 0
295     fVerV (i,j,1) = 0. _d 0
296     fVerV (i,j,2) = 0. _d 0
297 cnh 1.27 phiHyd(i,j,1) = 0. _d 0
298 cnh 1.30 K13 (i,j,1) = 0. _d 0
299     K23 (i,j,1) = 0. _d 0
300     K33 (i,j,1) = 0. _d 0
301     KapGM (i,j) = GMkbackground
302 cnh 1.7 ENDDO
303     ENDDO
304    
305 adcroft 1.45 DO k=1,Nr
306     DO j=1-OLy,sNy+OLy
307     DO i=1-OLx,sNx+OLx
308     #ifdef INCLUDE_CONVECT_CALL
309     ConvectCount(i,j,k) = 0.
310     #endif
311     KappaRT(i,j,k) = 0. _d 0
312     KappaRS(i,j,k) = 0. _d 0
313     ENDDO
314     ENDDO
315     ENDDO
316    
317 cnh 1.1 iMin = 1-OLx+1
318     iMax = sNx+OLx
319     jMin = 1-OLy+1
320     jMax = sNy+OLy
321 cnh 1.35
322 cnh 1.1
323 cnh 1.19 K = 1
324 cnh 1.31 BOTTOM_LAYER = K .EQ. Nr
325 cnh 1.19
326 cnh 1.38 #ifdef DO_PIPELINED_CORRECTION_STEP
327 adcroft 1.4 C-- Calculate gradient of surface pressure
328 cnh 1.28 CALL CALC_GRAD_ETA_SURF(
329 adcroft 1.4 I bi,bj,iMin,iMax,jMin,jMax,
330 cnh 1.29 O etaSurfX,etaSurfY,
331 adcroft 1.4 I myThid)
332     C-- Update fields in top level according to tendency terms
333 adcroft 1.11 CALL CORRECTION_STEP(
334 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K,
335     I etaSurfX,etaSurfY,myTime,myThid)
336 heimbach 1.49
337 adcroft 1.43 #ifdef ALLOW_OBCS
338 heimbach 1.49 IF (openBoundaries) THEN
339     #ifdef ALLOW_AUTODIFF_TAMC
340     CADJ STORE uvel (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
341     CADJ STORE vvel (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
342     CADJ STORE theta(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
343     CADJ STORE salt(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
344     #endif
345     CALL APPLY_OBCS1( bi, bj, K, myThid )
346     END IF
347 adcroft 1.43 #endif
348 heimbach 1.49
349 adcroft 1.21 IF ( .NOT. BOTTOM_LAYER ) THEN
350     C-- Update fields in layer below according to tendency terms
351     CALL CORRECTION_STEP(
352 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K+1,
353     I etaSurfX,etaSurfY,myTime,myThid)
354 adcroft 1.43 #ifdef ALLOW_OBCS
355 heimbach 1.49 IF (openBoundaries) THEN
356     #ifdef ALLOW_AUTODIFF_TAMC
357     CADJ STORE uvel (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
358     CADJ STORE vvel (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
359     CADJ STORE theta(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
360     CADJ STORE salt(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
361     #endif
362     CALL APPLY_OBCS1( bi, bj, K+1, myThid )
363     END IF
364 adcroft 1.43 #endif
365 adcroft 1.21 ENDIF
366 cnh 1.38 #endif
367 cnh 1.7 C-- Density of 1st level (below W(1)) reference to level 1
368 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
369 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
370     CADJ STORE theta(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
371     CADJ STORE salt (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
372     #endif
373 cnh 1.7 CALL FIND_RHO(
374 cnh 1.19 I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
375 cnh 1.7 O rhoKm1,
376     I myThid )
377 cnh 1.38 #endif
378 cnh 1.19
379 adcroft 1.42 IF ( (.NOT. BOTTOM_LAYER)
380     #ifdef ALLOW_KPP
381     & .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
382     #endif
383     & ) THEN
384 cnh 1.19 C-- Check static stability with layer below
385 cnh 1.30 C-- and mix as needed.
386 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
387 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
388     CADJ STORE theta(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
389     CADJ STORE salt (:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
390     #endif
391 cnh 1.19 CALL FIND_RHO(
392     I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
393     O rhoKp1,
394     I myThid )
395 cnh 1.38 #endif
396 heimbach 1.49
397 cnh 1.38 #ifdef INCLUDE_CONVECT_CALL
398 heimbach 1.49
399     #ifdef ALLOW_AUTODIFF_TAMC
400     CADJ STORE rhoKm1(:,:) = comlev1_2d, key = ikey, byte = isbyte
401     CADJ STORE rhoKp1(:,:) = comlev1_2d, key = ikey, byte = isbyte
402     #endif
403 cnh 1.19 CALL CONVECT(
404     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
405 adcroft 1.45 U ConvectCount,
406 cnh 1.19 I myTime,myIter,myThid)
407 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
408     CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)
409     CADJ & = comlev1_2d, key = ikey, byte = isbyte
410     CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)
411     CADJ & = comlev1_2d, key = ikey, byte = isbyte
412 cnh 1.38 #endif
413 heimbach 1.49
414     #endif
415    
416 adcroft 1.45 C-- Implicit Vertical Diffusion for Convection
417     IF (ivdc_kappa.NE.0.) CALL CALC_IVDC(
418     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
419     U ConvectCount, KappaRT, KappaRS,
420     I myTime,myIter,myThid)
421 heimbach 1.49 CRG: do we need do store STORE KappaRT, KappaRS ?
422    
423 cnh 1.19 C-- Recompute density after mixing
424 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
425 cnh 1.19 CALL FIND_RHO(
426     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
427     O rhoKm1,
428     I myThid )
429 cnh 1.38 #endif
430 cnh 1.19 ENDIF
431 cnh 1.26 C-- Calculate buoyancy
432 cnh 1.32 CALL CALC_BUOYANCY(
433 cnh 1.26 I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,
434     O buoyKm1,
435     I myThid )
436 cnh 1.38 C-- Integrate hydrostatic balance for phiHyd with BC of
437     C-- phiHyd(z=0)=0
438 cnh 1.26 CALL CALC_PHI_HYD(
439     I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1,
440     U phiHyd,
441 adcroft 1.5 I myThid )
442    
443 heimbach 1.49 C----------------------------------------------
444     C-- start of downward loop
445     C----------------------------------------------
446 cnh 1.31 DO K=2,Nr
447 heimbach 1.49
448     #ifdef ALLOW_AUTODIFF_TAMC
449     kkey = ikact*(Nr-2+1) + (k-2) + 1
450     #endif
451    
452 cnh 1.31 BOTTOM_LAYER = K .EQ. Nr
453 heimbach 1.49
454 cnh 1.38 #ifdef DO_PIPELINED_CORRECTION_STEP
455 adcroft 1.21 IF ( .NOT. BOTTOM_LAYER ) THEN
456     C-- Update fields in layer below according to tendency terms
457     CALL CORRECTION_STEP(
458 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K+1,
459     I etaSurfX,etaSurfY,myTime,myThid)
460 adcroft 1.43 #ifdef ALLOW_OBCS
461 heimbach 1.49 IF (openBoundaries) THEN
462     #ifdef ALLOW_AUTODIFF_TAMC
463     CADJ STORE uvel (:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
464     CADJ STORE vvel (:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
465     CADJ STORE theta(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
466     CADJ STORE salt(:,:,k,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
467     #endif
468     CALL APPLY_OBCS1( bi, bj, K+1, myThid )
469     END IF
470 adcroft 1.43 #endif
471 adcroft 1.21 ENDIF
472 cnh 1.38 #endif
473 heimbach 1.49
474 cnh 1.19 C-- Density of K level (below W(K)) reference to K level
475 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
476 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
477     CADJ STORE theta(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
478     CADJ STORE salt (:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
479     #endif
480 cnh 1.19 CALL FIND_RHO(
481     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
482     O rhoK,
483     I myThid )
484 cnh 1.38 #endif
485 adcroft 1.42 IF ( (.NOT. BOTTOM_LAYER)
486     #ifdef ALLOW_KPP
487     & .AND. (.NOT.usingKPPmixing) ! CONVECT not needed with KPP mixing
488     #endif
489     & ) THEN
490 cnh 1.27 C-- Check static stability with layer below and mix as needed.
491     C-- Density of K+1 level (below W(K+1)) reference to K level.
492 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
493 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
494     CADJ STORE theta(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
495     CADJ STORE salt (:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
496     #endif
497 cnh 1.19 CALL FIND_RHO(
498     I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType,
499     O rhoKp1,
500     I myThid )
501 cnh 1.38 #endif
502 heimbach 1.49
503     #ifdef ALLOW_AUTODIFF_TAMC
504     CADJ STORE rhok (:,:) = comlev1_3d, key = kkey, byte = isbyte
505     CADJ STORE rhoKm1(:,:) = comlev1_3d, key = kkey, byte = isbyte
506     CADJ STORE rhoKp1(:,:) = comlev1_3d, key = kkey, byte = isbyte
507     #endif
508    
509 cnh 1.38 #ifdef INCLUDE_CONVECT_CALL
510 cnh 1.19 CALL CONVECT(
511     I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1,
512 adcroft 1.45 U ConvectCount,
513 cnh 1.19 I myTime,myIter,myThid)
514 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
515     CADJ STORE theta(:,:,k+1,bi,bj),theta(:,:,k,bi,bj)
516     CADJ & = comlev1_3d, key = kkey, byte = isbyte
517     CADJ STORE salt (:,:,k+1,bi,bj),salt (:,:,k,bi,bj)
518     CADJ & = comlev1_3d, key = kkey, byte = isbyte
519     #endif
520 cnh 1.38 #endif
521 heimbach 1.49
522 adcroft 1.45 C-- Implicit Vertical Diffusion for Convection
523 heimbach 1.49 IF (ivdc_kappa.NE.0.) THEN
524     CALL CALC_IVDC(
525 adcroft 1.45 I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1,
526     U ConvectCount, KappaRT, KappaRS,
527     I myTime,myIter,myThid)
528 heimbach 1.49 CRG: do we need do store STORE KappaRT, KappaRS ?
529     END IF
530    
531 cnh 1.19 C-- Recompute density after mixing
532 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
533 cnh 1.19 CALL FIND_RHO(
534     I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType,
535     O rhoK,
536     I myThid )
537 cnh 1.38 #endif
538 cnh 1.19 ENDIF
539 cnh 1.26 C-- Calculate buoyancy
540 cnh 1.32 CALL CALC_BUOYANCY(
541 cnh 1.26 I bi,bj,iMin,iMax,jMin,jMax,K,rhoK,
542     O buoyK,
543     I myThid )
544 cnh 1.38 C-- Integrate hydrostatic balance for phiHyd with BC of
545     C-- phiHyd(z=0)=0
546 cnh 1.26 CALL CALC_PHI_HYD(
547 cnh 1.30 I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK,
548     U phiHyd,
549     I myThid )
550 cnh 1.19 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
551 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
552 cnh 1.19 CALL FIND_RHO(
553 cnh 1.30 I bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType,
554     O rhoTmp,
555     I myThid )
556 cnh 1.38 #endif
557 heimbach 1.49
558 cnh 1.38 #ifdef INCLUDE_CALC_ISOSLOPES_CALL
559 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
560     CADJ STORE rhoTmp(:,:) = comlev1_3d, key = kkey, byte = isbyte
561     CADJ STORE rhok (:,:) = comlev1_3d, key = kkey, byte = isbyte
562     CADJ STORE rhoKm1(:,:) = comlev1_3d, key = kkey, byte = isbyte
563     CADJ STORE kapgm (:,:) = comlev1_3d, key = kkey, byte = isbyte
564     #endif
565 cnh 1.19 CALL CALC_ISOSLOPES(
566 cnh 1.30 I bi, bj, iMin, iMax, jMin, jMax, K,
567     I rhoKm1, rhoK, rhotmp,
568     O K13, K23, K33, KapGM,
569     I myThid )
570 cnh 1.38 #endif
571 heimbach 1.49
572 cnh 1.19 DO J=jMin,jMax
573     DO I=iMin,iMax
574 cnh 1.38 #ifdef INCLUDE_FIND_RHO_CALL
575 cnh 1.27 rhoKm1 (I,J) = rhoK(I,J)
576 cnh 1.38 #endif
577 cnh 1.27 buoyKm1(I,J) = buoyK(I,J)
578 cnh 1.19 ENDDO
579 adcroft 1.10 ENDDO
580 heimbach 1.49 ENDDO
581     C-- end of k loop
582    
583     #ifdef ALLOW_AUTODIFF_TAMC
584     CADJ STORE theta(:,:,:,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
585     CADJ STORE salt (:,:,:,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
586     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
587     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
588     #endif
589 adcroft 1.11
590 adcroft 1.42 #ifdef ALLOW_KPP
591 heimbach 1.49 C----------------------------------------------
592 adcroft 1.42 C-- Compute KPP mixing coefficients
593 heimbach 1.49 C----------------------------------------------
594 adcroft 1.42 IF (usingKPPmixing) THEN
595 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
596     CADJ STORE fu (:,: ,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
597     CADJ STORE fv (:,: ,bi,bj) = comlev1_2d, key = ikey, byte = isbyte
598     #endif
599 adcroft 1.42 CALL TIMER_START('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
600     I , myThid)
601     CALL KVMIX(
602     I bi, bj, myTime, myThid )
603     CALL TIMER_STOP ('KVMIX (FIND KPP COEFFICIENTS) [DYNAMICS]'
604     I , myThid)
605     ENDIF
606     #endif
607    
608 heimbach 1.49 C----------------------------------------------
609     C-- start of upward loop
610     C----------------------------------------------
611 cnh 1.31 DO K = Nr, 1, -1
612 cnh 1.30
613 cnh 1.1 kM1 =max(1,k-1) ! Points to level above k (=k-1)
614     kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
615     kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer
616 heimbach 1.49
617 cnh 1.1 iMin = 1-OLx+2
618     iMax = sNx+OLx-1
619     jMin = 1-OLy+2
620     jMax = sNy+OLy-1
621    
622 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
623     kkey = ikact*(Nr-1+1) + (k-1) + 1
624     #endif
625    
626     #ifdef ALLOW_AUTODIFF_TAMC
627     CADJ STORE rvel (:,:,kDown) = comlev1_3d, key = kkey, byte = isbyte
628     CADJ STORE rTrans(:,:) = comlev1_3d, key = kkey, byte = isbyte
629     CADJ STORE KappaRT(:,:,:) = comlev1_3d, key = kkey, byte = isbyte
630     CADJ STORE KappaRS(:,:,:) = comlev1_3d, key = kkey, byte = isbyte
631     #endif
632    
633 cnh 1.1 C-- Get temporary terms used by tendency routines
634     CALL CALC_COMMON_FACTORS (
635     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
636 cnh 1.30 O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp,
637 cnh 1.1 I myThid)
638 heimbach 1.49
639 adcroft 1.47 #ifdef ALLOW_OBCS
640     IF (openBoundaries) THEN
641     CALL APPLY_OBCS3( bi, bj, K, Kup, rTrans, rVel, myThid )
642     ENDIF
643     #endif
644 heimbach 1.49
645 cnh 1.38 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
646 adcroft 1.12 C-- Calculate the total vertical diffusivity
647     CALL CALC_DIFFUSIVITY(
648     I bi,bj,iMin,iMax,jMin,jMax,K,
649     I maskC,maskUp,KapGM,K33,
650 adcroft 1.42 O KappaRT,KappaRS,KappaRU,KappaRV,
651 adcroft 1.12 I myThid)
652 cnh 1.38 #endif
653 cnh 1.1 C-- Calculate accelerations in the momentum equations
654 cnh 1.9 IF ( momStepping ) THEN
655     CALL CALC_MOM_RHS(
656     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
657 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,rVel,maskC,
658 adcroft 1.42 I phiHyd,KappaRU,KappaRV,
659 cnh 1.9 U aTerm,xTerm,cTerm,mTerm,pTerm,
660     U fZon, fMer, fVerU, fVerV,
661 cnh 1.38 I myTime, myThid)
662 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
663     #ifdef INCLUDE_CD_CODE
664     ELSE
665     DO j=1-OLy,sNy+OLy
666     DO i=1-OLx,sNx+OLx
667     guCD(i,j,k,bi,bj) = 0.0
668     gvCD(i,j,k,bi,bj) = 0.0
669     END DO
670     END DO
671     #endif
672     #endif
673 cnh 1.9 ENDIF
674 cnh 1.1 C-- Calculate active tracer tendencies
675 cnh 1.9 IF ( tempStepping ) THEN
676     CALL CALC_GT(
677     I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
678 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
679     I K13,K23,KappaRT,KapGM,
680 cnh 1.9 U aTerm,xTerm,fZon,fMer,fVerT,
681 cnh 1.37 I myTime, myThid)
682 cnh 1.9 ENDIF
683 adcroft 1.18 IF ( saltStepping ) THEN
684     CALL CALC_GS(
685     I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
686 cnh 1.30 I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC,
687     I K13,K23,KappaRS,KapGM,
688 adcroft 1.18 U aTerm,xTerm,fZon,fMer,fVerS,
689 cnh 1.37 I myTime, myThid)
690 adcroft 1.18 ENDIF
691 adcroft 1.47 #ifdef ALLOW_OBCS
692     C-- Calculate future values on open boundaries
693     IF (openBoundaries) THEN
694     Caja CALL CYCLE_OBCS( K, bi, bj, myThid )
695     CALL SET_OBCS( K, bi, bj, myTime+deltaTclock, myThid )
696     ENDIF
697     #endif
698 adcroft 1.11 C-- Prediction step (step forward all model variables)
699     CALL TIMESTEP(
700     I bi,bj,iMin,iMax,jMin,jMax,K,
701 adcroft 1.46 I myIter, myThid)
702 adcroft 1.43 #ifdef ALLOW_OBCS
703 adcroft 1.41 C-- Apply open boundary conditions
704 heimbach 1.49 IF (openBoundaries) THEN
705     #ifdef ALLOW_AUTODIFF_TAMC
706     CADJ STORE gunm1(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
707     CADJ STORE gvnm1(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
708     CADJ STORE gwnm1(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
709     #endif
710     CALL APPLY_OBCS2( bi, bj, K, myThid )
711     END IF
712 adcroft 1.43 #endif
713 adcroft 1.41 C-- Freeze water
714 heimbach 1.49 IF (allowFreezing) THEN
715     #ifdef ALLOW_AUTODIFF_TAMC
716     CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_3d, key = kkey, byte = isbyte
717     #endif
718     CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, K, myThid )
719     END IF
720 adcroft 1.48
721     #ifdef DIVG_IN_DYNAMICS
722 adcroft 1.11 C-- Diagnose barotropic divergence of predicted fields
723 cnh 1.31 CALL CALC_DIV_GHAT(
724 adcroft 1.11 I bi,bj,iMin,iMax,jMin,jMax,K,
725     I xA,yA,
726     I myThid)
727 adcroft 1.48 #endif /* DIVG_IN_DYNAMICS */
728 adcroft 1.23
729     C-- Cumulative diagnostic calculations (ie. time-averaging)
730 cnh 1.38 #ifdef INCLUDE_DIAGNOSTICS_INTERFACE_CODE
731 adcroft 1.23 IF (taveFreq.GT.0.) THEN
732     CALL DO_TIME_AVERAGES(
733     I myTime, myIter, bi, bj, K, kUp, kDown,
734 adcroft 1.45 I K13, K23, rVel, KapGM, ConvectCount,
735 adcroft 1.23 I myThid )
736     ENDIF
737     #endif
738 adcroft 1.45
739 adcroft 1.11
740     ENDDO ! K
741 adcroft 1.12
742     C-- Implicit diffusion
743     IF (implicitDiffusion) THEN
744 heimbach 1.49
745     #ifdef ALLOW_AUTODIFF_TAMC
746     maximpl = 6
747     iikey = ikact*maximpl
748     #endif
749    
750     IF (tempStepping) THEN
751     #ifdef ALLOW_AUTODIFF_TAMC
752     idkey = iikey + 1
753     #endif
754     CALL IMPLDIFF(
755 adcroft 1.42 I bi, bj, iMin, iMax, jMin, jMax,
756     I deltaTtracer, KappaRT,recip_HFacC,
757     U gTNm1,
758     I myThid )
759 heimbach 1.49 END IF
760    
761     IF (saltStepping) THEN
762     #ifdef ALLOW_AUTODIFF_TAMC
763     idkey = iikey + 2
764     #endif
765     CALL IMPLDIFF(
766 adcroft 1.42 I bi, bj, iMin, iMax, jMin, jMax,
767     I deltaTtracer, KappaRS,recip_HFacC,
768     U gSNm1,
769     I myThid )
770 heimbach 1.49 END IF
771    
772 adcroft 1.44 ENDIF ! implicitDiffusion
773 heimbach 1.49
774 adcroft 1.44 C-- Implicit viscosity
775     IF (implicitViscosity) THEN
776 heimbach 1.49
777 adcroft 1.42 IF (momStepping) THEN
778 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
779     idkey = iikey + 3
780     #endif
781 adcroft 1.42 CALL IMPLDIFF(
782     I bi, bj, iMin, iMax, jMin, jMax,
783     I deltaTmom, KappaRU,recip_HFacW,
784     U gUNm1,
785     I myThid )
786 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
787     idkey = iikey + 4
788     #endif
789 adcroft 1.42 CALL IMPLDIFF(
790     I bi, bj, iMin, iMax, jMin, jMax,
791     I deltaTmom, KappaRV,recip_HFacS,
792     U gVNm1,
793     I myThid )
794 heimbach 1.49
795 adcroft 1.42 #ifdef INCLUDE_CD_CODE
796 heimbach 1.49
797     #ifdef ALLOW_AUTODIFF_TAMC
798     idkey = iikey + 5
799     #endif
800 adcroft 1.42 CALL IMPLDIFF(
801     I bi, bj, iMin, iMax, jMin, jMax,
802     I deltaTmom, KappaRU,recip_HFacW,
803     U vVelD,
804     I myThid )
805 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
806     idkey = iikey + 6
807     #endif
808 adcroft 1.42 CALL IMPLDIFF(
809     I bi, bj, iMin, iMax, jMin, jMax,
810     I deltaTmom, KappaRV,recip_HFacS,
811     U uVelD,
812     I myThid )
813 heimbach 1.49
814 adcroft 1.42 #endif
815 heimbach 1.49
816 adcroft 1.42 ENDIF ! momStepping
817 adcroft 1.44 ENDIF ! implicitViscosity
818 cnh 1.1
819     ENDDO
820     ENDDO
821 adcroft 1.6
822 cnh 1.19 C write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)),
823     C & maxval(cg2d_x(1:sNx,1:sNy,:,:))
824 cnh 1.20 C write(0,*) 'dynamics: U ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.),
825 adcroft 1.21 C & maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.)
826 cnh 1.20 C write(0,*) 'dynamics: V ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.),
827 adcroft 1.21 C & maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.)
828 cnh 1.30 C write(0,*) 'dynamics: rVel(1) ',
829     C & minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.),
830     C & maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.)
831     C write(0,*) 'dynamics: rVel(2) ',
832     C & minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.),
833     C & maxval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.)
834 adcroft 1.15 cblk write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
835     cblk & maxval(K13(1:sNx,1:sNy,:))
836     cblk write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
837     cblk & maxval(K23(1:sNx,1:sNy,:))
838     cblk write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
839     cblk & maxval(K33(1:sNx,1:sNy,:))
840 cnh 1.19 C write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)),
841     C & maxval(gT(1:sNx,1:sNy,:,:,:))
842     C write(0,*) 'dynamics: T ',minval(Theta(1:sNx,1:sNy,:,:,:)),
843     C & maxval(Theta(1:sNx,1:sNy,:,:,:))
844     C write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)),
845     C & maxval(gS(1:sNx,1:sNy,:,:,:))
846     C write(0,*) 'dynamics: S ',minval(salt(1:sNx,1:sNy,:,:,:)),
847     C & maxval(salt(1:sNx,1:sNy,:,:,:))
848 cnh 1.31 C write(0,*) 'dynamics: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.),
849     C & maxval(phiHyd/(Gravity*Rhonil))
850 cnh 1.36 C CALL PLOT_FIELD_XYZRL( gU, ' GU exiting dyanmics ' ,
851     C &Nr, 1, myThid )
852     C CALL PLOT_FIELD_XYZRL( gV, ' GV exiting dyanmics ' ,
853     C &Nr, 1, myThid )
854     C CALL PLOT_FIELD_XYZRL( gS, ' GS exiting dyanmics ' ,
855     C &Nr, 1, myThid )
856     C CALL PLOT_FIELD_XYZRL( gT, ' GT exiting dyanmics ' ,
857 cnh 1.38 C &Nr, 1, myThid )
858     C CALL PLOT_FIELD_XYZRL( phiHyd, ' phiHyd exiting dyanmics ' ,
859 cnh 1.36 C &Nr, 1, myThid )
860    
861 cnh 1.1
862     RETURN
863     END

  ViewVC Help
Powered by ViewVC 1.1.22