/[MITgcm]/MITgcm_contrib/AITCZ/code/dynamics.F
ViewVC logotype

Annotation of /MITgcm_contrib/AITCZ/code/dynamics.F

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


Revision 1.1 - (hide annotations) (download)
Wed Aug 20 15:24:59 2003 UTC (21 years, 11 months ago) by czaja
Branch: MAIN
CVS Tags: HEAD
Initial creation of Arnaud's simple coupled simulation.

1 czaja 1.1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/dynamics.F,v 1.83 2001/09/27 20:12:10 heimbach Exp $
2     C $Name: release1_beta1 $
3    
4     #include "CPP_OPTIONS.h"
5    
6     CBOP
7     C !ROUTINE: DYNAMICS
8     C !INTERFACE:
9     SUBROUTINE DYNAMICS(myTime, myIter, myThid)
10     C !DESCRIPTION: \bv
11     C *==========================================================*
12     C | SUBROUTINE DYNAMICS
13     C | o Controlling routine for the explicit part of the model
14     C | dynamics.
15     C *==========================================================*
16     C | This routine evaluates the "dynamics" terms for each
17     C | block of ocean in turn. Because the blocks of ocean have
18     C | overlap regions they are independent of one another.
19     C | If terms involving lateral integrals are needed in this
20     C | routine care will be needed. Similarly finite-difference
21     C | operations with stencils wider than the overlap region
22     C | require special consideration.
23     C | The algorithm...
24     C |
25     C | "Correction Step"
26     C | =================
27     C | Here we update the horizontal velocities with the surface
28     C | pressure such that the resulting flow is either consistent
29     C | with the free-surface evolution or the rigid-lid:
30     C | U[n] = U* + dt x d/dx P
31     C | V[n] = V* + dt x d/dy P
32     C |
33     C | "Calculation of Gs"
34     C | ===================
35     C | This is where all the accelerations and tendencies (ie.
36     C | physics, parameterizations etc...) are calculated
37     C | rho = rho ( theta[n], salt[n] )
38     C | b = b(rho, theta)
39     C | K31 = K31 ( rho )
40     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
41     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
42     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
43     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
44     C |
45     C | "Time-stepping" or "Prediction"
46     C | ================================
47     C | The models variables are stepped forward with the appropriate
48     C | time-stepping scheme (currently we use Adams-Bashforth II)
49     C | - For momentum, the result is always *only* a "prediction"
50     C | in that the flow may be divergent and will be "corrected"
51     C | later with a surface pressure gradient.
52     C | - Normally for tracers the result is the new field at time
53     C | level [n+1} *BUT* in the case of implicit diffusion the result
54     C | is also *only* a prediction.
55     C | - We denote "predictors" with an asterisk (*).
56     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
57     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
58     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
59     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
60     C | With implicit diffusion:
61     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63     C | (1 + dt * K * d_zz) theta[n] = theta*
64     C | (1 + dt * K * d_zz) salt[n] = salt*
65     C |
66     C *==========================================================*
67     C \ev
68     C !USES:
69     IMPLICIT NONE
70     C == Global variables ===
71     #include "SIZE.h"
72     #include "EEPARAMS.h"
73     #include "PARAMS.h"
74     #include "DYNVARS.h"
75     #include "GRID.h"
76     #ifdef ALLOW_PASSIVE_TRACER
77     #include "TR1.h"
78     #endif
79     #ifdef ALLOW_AUTODIFF_TAMC
80     # include "tamc.h"
81     # include "tamc_keys.h"
82     # include "FFIELDS.h"
83     # ifdef ALLOW_KPP
84     # include "KPP.h"
85     # endif
86     # ifdef ALLOW_GMREDI
87     # include "GMREDI.h"
88     # endif
89     #endif /* ALLOW_AUTODIFF_TAMC */
90     #ifdef ALLOW_TIMEAVE
91     #include "TIMEAVE_STATV.h"
92     #endif
93    
94     C !CALLING SEQUENCE:
95     C DYNAMICS()
96     C |
97     C |-- CALC_GRAD_PHI_SURF
98     C |
99     C |-- CALC_VISCOSITY
100     C |
101     C |-- CALC_PHI_HYD
102     C |
103     C |-- MOM_FLUXFORM
104     C |
105     C |-- MOM_VECINV
106     C |
107     C |-- TIMESTEP
108     C |
109     C |-- OBCS_APPLY_UV
110     C |
111     C |-- IMPLDIFF
112     C |
113     C |-- OBCS_APPLY_UV
114     C |
115     C |-- CALL TIMEAVE_CUMUL_1T
116     C |-- CALL TIMEAVE_CUMULATE
117     C |-- CALL DEBUG_STATS_RL
118    
119     C !INPUT/OUTPUT PARAMETERS:
120     C == Routine arguments ==
121     C myTime - Current time in simulation
122     C myIter - Current iteration number in simulation
123     C myThid - Thread number for this instance of the routine.
124     _RL myTime
125     INTEGER myIter
126     INTEGER myThid
127    
128     C !LOCAL VARIABLES:
129     C == Local variables
130     C fVer[STUV] o fVer: Vertical flux term - note fVer
131     C is "pipelined" in the vertical
132     C so we need an fVer for each
133     C variable.
134     C rhoK, rhoKM1 - Density at current level, and level above
135     C phiHyd - Hydrostatic part of the potential phiHydi.
136     C In z coords phiHydiHyd is the hydrostatic
137     C Potential (=pressure/rho0) anomaly
138     C In p coords phiHydiHyd is the geopotential
139     C surface height anomaly.
140     C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
141     C phiSurfY or geopotentiel (atmos) in X and Y direction
142     C iMin, iMax - Ranges and sub-block indices on which calculations
143     C jMin, jMax are applied.
144     C bi, bj
145     C k, kup, - Index for layer above and below. kup and kDown
146     C kDown, km1 are switched with layer to be the appropriate
147     C index into fVerTerm.
148     _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
150     _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
151     _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153     _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155     _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
157     _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
158     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
159     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
160    
161     C This is currently used by IVDC and Diagnostics
162     _RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
163    
164     INTEGER iMin, iMax
165     INTEGER jMin, jMax
166     INTEGER bi, bj
167     INTEGER i, j
168     INTEGER k, km1, kp1, kup, kDown
169    
170     Cjmc : add for phiHyd output <- but not working if multi tile per CPU
171     c CHARACTER*(MAX_LEN_MBUF) suff
172     c LOGICAL DIFFERENT_MULTIPLE
173     c EXTERNAL DIFFERENT_MULTIPLE
174     Cjmc(end)
175    
176     C--- The algorithm...
177     C
178     C "Correction Step"
179     C =================
180     C Here we update the horizontal velocities with the surface
181     C pressure such that the resulting flow is either consistent
182     C with the free-surface evolution or the rigid-lid:
183     C U[n] = U* + dt x d/dx P
184     C V[n] = V* + dt x d/dy P
185     C
186     C "Calculation of Gs"
187     C ===================
188     C This is where all the accelerations and tendencies (ie.
189     C physics, parameterizations etc...) are calculated
190     C rho = rho ( theta[n], salt[n] )
191     C b = b(rho, theta)
192     C K31 = K31 ( rho )
193     C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
194     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
195     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
196     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
197     C
198     C "Time-stepping" or "Prediction"
199     C ================================
200     C The models variables are stepped forward with the appropriate
201     C time-stepping scheme (currently we use Adams-Bashforth II)
202     C - For momentum, the result is always *only* a "prediction"
203     C in that the flow may be divergent and will be "corrected"
204     C later with a surface pressure gradient.
205     C - Normally for tracers the result is the new field at time
206     C level [n+1} *BUT* in the case of implicit diffusion the result
207     C is also *only* a prediction.
208     C - We denote "predictors" with an asterisk (*).
209     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
210     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
211     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
212     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213     C With implicit diffusion:
214     C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
215     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
216     C (1 + dt * K * d_zz) theta[n] = theta*
217     C (1 + dt * K * d_zz) salt[n] = salt*
218     C---
219     CEOP
220    
221     C-- Set up work arrays with valid (i.e. not NaN) values
222     C These inital values do not alter the numerical results. They
223     C just ensure that all memory references are to valid floating
224     C point numbers. This prevents spurious hardware signals due to
225     C uninitialised but inert locations.
226     DO j=1-OLy,sNy+OLy
227     DO i=1-OLx,sNx+OLx
228     DO k=1,Nr
229     phiHyd(i,j,k) = 0. _d 0
230     KappaRU(i,j,k) = 0. _d 0
231     KappaRV(i,j,k) = 0. _d 0
232     sigmaX(i,j,k) = 0. _d 0
233     sigmaY(i,j,k) = 0. _d 0
234     sigmaR(i,j,k) = 0. _d 0
235     ENDDO
236     rhoKM1 (i,j) = 0. _d 0
237     rhok (i,j) = 0. _d 0
238     phiSurfX(i,j) = 0. _d 0
239     phiSurfY(i,j) = 0. _d 0
240     ENDDO
241     ENDDO
242    
243     #ifdef ALLOW_AUTODIFF_TAMC
244     C-- HPF directive to help TAMC
245     CHPF$ INDEPENDENT
246     #endif /* ALLOW_AUTODIFF_TAMC */
247    
248     DO bj=myByLo(myThid),myByHi(myThid)
249    
250     #ifdef ALLOW_AUTODIFF_TAMC
251     C-- HPF directive to help TAMC
252     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
253     CHPF$& ,phiHyd
254     CHPF$& ,KappaRU,KappaRV
255     CHPF$& )
256     #endif /* ALLOW_AUTODIFF_TAMC */
257    
258     DO bi=myBxLo(myThid),myBxHi(myThid)
259    
260     #ifdef ALLOW_AUTODIFF_TAMC
261     act1 = bi - myBxLo(myThid)
262     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
263     act2 = bj - myByLo(myThid)
264     max2 = myByHi(myThid) - myByLo(myThid) + 1
265     act3 = myThid - 1
266     max3 = nTx*nTy
267     act4 = ikey_dynamics - 1
268     ikey = (act1 + 1) + act2*max1
269     & + act3*max1*max2
270     & + act4*max1*max2*max3
271     #endif /* ALLOW_AUTODIFF_TAMC */
272    
273     C-- Set up work arrays that need valid initial values
274     DO j=1-OLy,sNy+OLy
275     DO i=1-OLx,sNx+OLx
276     fVerU (i,j,1) = 0. _d 0
277     fVerU (i,j,2) = 0. _d 0
278     fVerV (i,j,1) = 0. _d 0
279     fVerV (i,j,2) = 0. _d 0
280     ENDDO
281     ENDDO
282    
283     C-- Start computation of dynamics
284     iMin = 1-OLx+2
285     iMax = sNx+OLx-1
286     jMin = 1-OLy+2
287     jMax = sNy+OLy-1
288    
289     #ifdef ALLOW_AUTODIFF_TAMC
290     CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte
291     #endif /* ALLOW_AUTODIFF_TAMC */
292    
293     C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
294     C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
295     IF (implicSurfPress.NE.1.) THEN
296     CALL CALC_GRAD_PHI_SURF(
297     I bi,bj,iMin,iMax,jMin,jMax,
298     I etaN,
299     O phiSurfX,phiSurfY,
300     I myThid )
301     ENDIF
302    
303     #ifdef ALLOW_AUTODIFF_TAMC
304     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
305     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
306     #ifdef ALLOW_KPP
307     CADJ STORE KPPviscAz (:,:,:,bi,bj)
308     CADJ & = comlev1_bibj, key=ikey, byte=isbyte
309     #endif /* ALLOW_KPP */
310     #endif /* ALLOW_AUTODIFF_TAMC */
311    
312     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
313     C-- Calculate the total vertical diffusivity
314     DO k=1,Nr
315     CALL CALC_VISCOSITY(
316     I bi,bj,iMin,iMax,jMin,jMax,k,
317     O KappaRU,KappaRV,
318     I myThid)
319     ENDDO
320     #endif
321    
322     C-- Start of dynamics loop
323     DO k=1,Nr
324    
325     C-- km1 Points to level above k (=k-1)
326     C-- kup Cycles through 1,2 to point to layer above
327     C-- kDown Cycles through 2,1 to point to current layer
328    
329     km1 = MAX(1,k-1)
330     kp1 = MIN(k+1,Nr)
331     kup = 1+MOD(k+1,2)
332     kDown= 1+MOD(k,2)
333    
334     #ifdef ALLOW_AUTODIFF_TAMC
335     kkey = (ikey-1)*Nr + k
336     #endif /* ALLOW_AUTODIFF_TAMC */
337    
338     C-- Integrate hydrostatic balance for phiHyd with BC of
339     C phiHyd(z=0)=0
340     C distinguishe between Stagger and Non Stagger time stepping
341     IF (staggerTimeStep) THEN
342     CALL CALC_PHI_HYD(
343     I bi,bj,iMin,iMax,jMin,jMax,k,
344     I gT, gS,
345     U phiHyd,
346     I myThid )
347     ELSE
348     CALL CALC_PHI_HYD(
349     I bi,bj,iMin,iMax,jMin,jMax,k,
350     I theta, salt,
351     U phiHyd,
352     I myThid )
353     ENDIF
354    
355     C-- Calculate accelerations in the momentum equations (gU, gV, ...)
356     C and step forward storing the result in gUnm1, gVnm1, etc...
357     IF ( momStepping ) THEN
358     #ifndef DISABLE_MOM_FLUXFORM
359     IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
360     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
361     I phiHyd,KappaRU,KappaRV,
362     U fVerU, fVerV,
363     I myTime, myIter, myThid)
364     #endif
365     #ifndef DISABLE_MOM_VECINV
366     IF (vectorInvariantMomentum) CALL MOM_VECINV(
367     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
368     I phiHyd,KappaRU,KappaRV,
369     U fVerU, fVerV,
370     I myTime, myIter, myThid)
371     #endif
372     CALL TIMESTEP(
373     I bi,bj,iMin,iMax,jMin,jMax,k,
374     I phiHyd, phiSurfX, phiSurfY,
375     I myIter, myThid)
376    
377     #ifdef ALLOW_OBCS
378     C-- Apply open boundary conditions
379     IF (useOBCS) THEN
380     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
381     END IF
382     #endif /* ALLOW_OBCS */
383    
384     #ifdef ALLOW_AUTODIFF_TAMC
385     #ifdef INCLUDE_CD_CODE
386     ELSE
387     DO j=1-OLy,sNy+OLy
388     DO i=1-OLx,sNx+OLx
389     guCD(i,j,k,bi,bj) = 0.0
390     gvCD(i,j,k,bi,bj) = 0.0
391     END DO
392     END DO
393     #endif /* INCLUDE_CD_CODE */
394     #endif /* ALLOW_AUTODIFF_TAMC */
395     ENDIF
396    
397    
398     C-- end of dynamics k loop (1:Nr)
399     ENDDO
400    
401    
402    
403     C-- Implicit viscosity
404     IF (implicitViscosity.AND.momStepping) THEN
405     #ifdef ALLOW_AUTODIFF_TAMC
406     idkey = iikey + 3
407     CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
408     #endif /* ALLOW_AUTODIFF_TAMC */
409     CALL IMPLDIFF(
410     I bi, bj, iMin, iMax, jMin, jMax,
411     I deltaTmom, KappaRU,recip_HFacW,
412     U gUNm1,
413     I myThid )
414     #ifdef ALLOW_AUTODIFF_TAMC
415     idkey = iikey + 4
416     CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
417     #endif /* ALLOW_AUTODIFF_TAMC */
418     CALL IMPLDIFF(
419     I bi, bj, iMin, iMax, jMin, jMax,
420     I deltaTmom, KappaRV,recip_HFacS,
421     U gVNm1,
422     I myThid )
423    
424     #ifdef ALLOW_OBCS
425     C-- Apply open boundary conditions
426     IF (useOBCS) THEN
427     DO K=1,Nr
428     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
429     ENDDO
430     END IF
431     #endif /* ALLOW_OBCS */
432    
433     #ifdef INCLUDE_CD_CODE
434     #ifdef ALLOW_AUTODIFF_TAMC
435     idkey = iikey + 5
436     CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
437     #endif /* ALLOW_AUTODIFF_TAMC */
438     CALL IMPLDIFF(
439     I bi, bj, iMin, iMax, jMin, jMax,
440     I deltaTmom, KappaRU,recip_HFacW,
441     U vVelD,
442     I myThid )
443     #ifdef ALLOW_AUTODIFF_TAMC
444     idkey = iikey + 6
445     CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte
446     #endif /* ALLOW_AUTODIFF_TAMC */
447     CALL IMPLDIFF(
448     I bi, bj, iMin, iMax, jMin, jMax,
449     I deltaTmom, KappaRV,recip_HFacS,
450     U uVelD,
451     I myThid )
452     #endif /* INCLUDE_CD_CODE */
453     C-- End If implicitViscosity.AND.momStepping
454     ENDIF
455    
456     Cjmc : add for phiHyd output <- but not working if multi tile per CPU
457     c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime)
458     c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN
459     c WRITE(suff,'(I10.10)') myIter+1
460     c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid)
461     c ENDIF
462     Cjmc(end)
463    
464     #ifdef ALLOW_TIMEAVE
465     IF (taveFreq.GT.0.) THEN
466     CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
467     I deltaTclock, bi, bj, myThid)
468     IF (ivdc_kappa.NE.0.) THEN
469     CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr,
470     I deltaTclock, bi, bj, myThid)
471     ENDIF
472     ENDIF
473     #endif /* ALLOW_TIMEAVE */
474    
475     ENDDO
476     ENDDO
477    
478     #ifndef DISABLE_DEBUGMODE
479     If (debugMode) THEN
480     CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
481     CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
482     CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
483     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
484     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
485     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
486     CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
487     CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
488     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
489     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
490     CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
491     CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
492     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
493     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
494     ENDIF
495     #endif
496    
497     RETURN
498     END

  ViewVC Help
Powered by ViewVC 1.1.22