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

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

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

revision 1.6 by adcroft, Wed May 20 21:29:31 1998 UTC revision 1.13 by adcroft, Mon Jun 8 18:45:28 1998 UTC
# Line 2  C $Header$ Line 2  C $Header$
2    
3  #include "CPP_EEOPTIONS.h"  #include "CPP_EEOPTIONS.h"
4    
5        SUBROUTINE DYNAMICS(myThid)        SUBROUTINE DYNAMICS(myTime, myIter, myThid)
6  C     /==========================================================\  C     /==========================================================\
7  C     | SUBROUTINE DYNAMICS                                      |  C     | SUBROUTINE DYNAMICS                                      |
8  C     | o Controlling routine for the explicit part of the model |  C     | o Controlling routine for the explicit part of the model |
# Line 29  C     == Global variables === Line 29  C     == Global variables ===
29  #include "DYNVARS.h"  #include "DYNVARS.h"
30    
31  C     == Routine arguments ==  C     == Routine arguments ==
32    C     myTime - Current time in simulation
33    C     myIter - Current iteration number in simulation
34  C     myThid - Thread number for this instance of the routine.  C     myThid - Thread number for this instance of the routine.
35        INTEGER myThid        INTEGER myThid
36          _RL myTime
37          INTEGER myIter
38    
39  C     == Local variables  C     == Local variables
40  C     xA, yA                 - Per block temporaries holding face areas  C     xA, yA                 - Per block temporaries holding face areas
# Line 80  C                          into fVerTerm Line 84  C                          into fVerTerm
84        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL pH    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
85        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87          _RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
91        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
92        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL K33   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)
93        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94          _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
95    
96        INTEGER iMin, iMax        INTEGER iMin, iMax
97        INTEGER jMin, jMax        INTEGER jMin, jMax
98        INTEGER bi, bj        INTEGER bi, bj
99        INTEGER i, j        INTEGER i, j
100        INTEGER k, kM1, kUp, kDown        INTEGER k, kM1, kUp, kDown
101    
102    C---    The algorithm...
103    C
104    C       "Correction Step"
105    C       =================
106    C       Here we update the horizontal velocities with the surface
107    C       pressure such that the resulting flow is either consistent
108    C       with the free-surface evolution or the rigid-lid:
109    C         U[n] = U* + dt x d/dx P
110    C         V[n] = V* + dt x d/dy P
111    C
112    C       "Calculation of Gs"
113    C       ===================
114    C       This is where all the accelerations and tendencies (ie.
115    C       physics, parameterizations etc...) are calculated
116    C         w = sum_z ( div. u[n] )
117    C         rho = rho ( theta[n], salt[n] )
118    C         K31 = K31 ( rho )
119    C         Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... )
120    C         Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... )
121    C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )
122    C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )
123    C
124    C       "Time-stepping" or "Prediction"
125    C       ================================
126    C       The models variables are stepped forward with the appropriate
127    C       time-stepping scheme (currently we use Adams-Bashforth II)
128    C       - For momentum, the result is always *only* a "prediction"
129    C       in that the flow may be divergent and will be "corrected"
130    C       later with a surface pressure gradient.
131    C       - Normally for tracers the result is the new field at time
132    C       level [n+1} *BUT* in the case of implicit diffusion the result
133    C       is also *only* a prediction.
134    C       - We denote "predictors" with an asterisk (*).
135    C         U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
136    C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
137    C         theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
138    C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
139    C       With implicit diffusion:
140    C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
141    C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
142    C         (1 + dt * K * d_zz) theta[n] = theta*
143    C         (1 + dt * K * d_zz) salt[n] = salt*
144    C---
145    
146  C--   Set up work arrays with valid (i.e. not NaN) values  C--   Set up work arrays with valid (i.e. not NaN) values
147  C     These inital values do not alter the numerical results. They  C     These inital values do not alter the numerical results. They
148  C     just ensure that all memory references are to valid floating  C     just ensure that all memory references are to valid floating
# Line 115  C     uninitialised but inert locations. Line 166  C     uninitialised but inert locations.
166           K13(i,j,k) = 0. _d 0           K13(i,j,k) = 0. _d 0
167           K23(i,j,k) = 0. _d 0           K23(i,j,k) = 0. _d 0
168           K33(i,j,k) = 0. _d 0           K33(i,j,k) = 0. _d 0
169             KappaZT(i,j,k) = 0. _d 0
170          ENDDO          ENDDO
171          rhokm1(i,j)  = 0. _d 0          rhokm1(i,j)  = 0. _d 0
172          rhokp1(i,j)  = 0. _d 0          rhokp1(i,j)  = 0. _d 0
173         ENDDO          rhotmp(i,j)  = 0. _d 0
       ENDDO  
 C--   Set up work arrays that need valid initial values  
       DO j=1-OLy,sNy+OLy  
        DO i=1-OLx,sNx+OLx  
         wTrans(i,j)  = 0. _d 0  
         fVerT(i,j,1) = 0. _d 0  
         fVerT(i,j,2) = 0. _d 0  
         fVerS(i,j,1) = 0. _d 0  
         fVerS(i,j,2) = 0. _d 0  
         fVerU(i,j,1) = 0. _d 0  
         fVerU(i,j,2) = 0. _d 0  
         fVerV(i,j,1) = 0. _d 0  
         fVerV(i,j,2) = 0. _d 0  
174         ENDDO         ENDDO
175        ENDDO        ENDDO
176    
177        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
178         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
179    
180  C--   Boundary condition on hydrostatic pressure is pH(z=0)=0  C--     Set up work arrays that need valid initial values
181          DO j=1-OLy,sNy+OLy          DO j=1-OLy,sNy+OLy
182           DO i=1-OLx,sNx+OLx           DO i=1-OLx,sNx+OLx
183              wTrans(i,j)  = 0. _d 0
184              fVerT(i,j,1) = 0. _d 0
185              fVerT(i,j,2) = 0. _d 0
186              fVerS(i,j,1) = 0. _d 0
187              fVerS(i,j,2) = 0. _d 0
188              fVerU(i,j,1) = 0. _d 0
189              fVerU(i,j,2) = 0. _d 0
190              fVerV(i,j,1) = 0. _d 0
191              fVerV(i,j,2) = 0. _d 0
192            pH(i,j,1) = 0. _d 0            pH(i,j,1) = 0. _d 0
193            K13(i,j,1) = 0. _d 0            K13(i,j,1) = 0. _d 0
194            K23(i,j,1) = 0. _d 0            K23(i,j,1) = 0. _d 0
# Line 161  C--     Calculate gradient of surface pr Line 209  C--     Calculate gradient of surface pr
209       I       myThid)       I       myThid)
210    
211  C--     Update fields in top level according to tendency terms  C--     Update fields in top level according to tendency terms
212          CALL TIMESTEP(          CALL CORRECTION_STEP(
213       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,1,pSurfX,pSurfY,myThid)
214    
215  C Density of 1st level (below W(1)) reference to level 1  C--     Density of 1st level (below W(1)) reference to level 1
216           CALL FIND_RHO(          CALL FIND_RHO(
217       I      bi, bj, iMin, iMax, jMin, jMax, 1, 1, 'LINEAR',       I     bi, bj, iMin, iMax, jMin, jMax, 1, 1, eosType,
218       O      rhoKm1,       O     rhoKm1,
219       I      myThid )       I     myThid )
220  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0
221           CALL CALC_PH(          CALL CALC_PH(
222       I       bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,       I      bi,bj,iMin,iMax,jMin,jMax,1,rhoKm1,rhoKm1,
223       U       pH,       U      pH,
224       I       myThid )       I      myThid )
225            DO J=1-Oly,sNy+Oly
226             DO I=1-Olx,sNx+Olx
227              rhoKp1(I,J)=rhoKm1(I,J)
228             ENDDO
229            ENDDO
230    
231          DO K=2,Nz          DO K=2,Nz
232  C--     Update fields in Kth level according to tendency terms  C--     Update fields in Kth level according to tendency terms
233          CALL TIMESTEP(          CALL CORRECTION_STEP(
234       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)       I       bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myThid)
235  C Density of K-1 level (above W(K)) reference to K level  C--     Density of K-1 level (above W(K)) reference to K-1 level
236           CALL FIND_RHO(  copt    CALL FIND_RHO(
237       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K, 'LINEAR',  copt I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,
238       O      rhoKm1,  copt O     rhoKm1,
239       I      myThid )  copt I     myThid )
240  C Density of K level (below W(K)) reference to K level  C       rhoKm1=rhoKp1
241           CALL FIND_RHO(          DO J=1-Oly,sNy+Oly
242       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',           DO I=1-Olx,sNx+Olx
243       O      rhoKp1,            rhoKm1(I,J)=rhoKp1(I,J)
244       I      myThid )           ENDDO
245            ENDDO
246    C--     Density of K level (below W(K)) reference to K level
247            CALL FIND_RHO(
248         I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
249         O     rhoKp1,
250         I     myThid )
251    C--     Density of K-1 level (above W(K)) reference to K level
252            CALL FIND_RHO(
253         I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K, eosType,
254         O     rhotmp,
255         I     myThid )
256  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation  C--     Calculate iso-neutral slopes for the GM/Redi parameterisation
257           CALL CALC_ISOSLOPES(          CALL CALC_ISOSLOPES(
258       I             bi, bj, iMin, iMax, jMin, jMax, K,       I            bi, bj, iMin, iMax, jMin, jMax, K,
259       I             rhoKm1, rhoKp1,       I            rhoKm1, rhoKp1, rhotmp,
260       O             K13, K23, K33, KapGM,       O            K13, K23, K33, KapGM,
261       I             myThid )       I            myThid )
262  C--     Calculate static stability and mix where convectively unstable  C--     Calculate static stability and mix where convectively unstable
263           CALL CONVECT(          CALL CONVECT(
264       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,myThid)       I      bi,bj,iMin,iMax,jMin,jMax,K,rhotmp,rhoKp1,
265  C Density of K-1 level (above W(K)) reference to K-1 level       I      myTime,myIter,myThid)
266           CALL FIND_RHO(  C--     Density of K-1 level (above W(K)) reference to K-1 level
267       I      bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, 'LINEAR',          CALL FIND_RHO(
268       O      rhoKm1,       I     bi, bj, iMin, iMax, jMin, jMax,  K-1, K-1, eosType,
269       I      myThid )       O     rhoKm1,
270  C Density of K level (below W(K)) referenced to K level       I     myThid )
271           CALL FIND_RHO(  C--     Density of K level (below W(K)) referenced to K level
272       I      bi, bj, iMin, iMax, jMin, jMax,  K, K, 'LINEAR',          CALL FIND_RHO(
273       O      rhoKp1,       I     bi, bj, iMin, iMax, jMin, jMax,  K, K, eosType,
274       I      myThid )       O     rhoKp1,
275         I     myThid )
276  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0  C--     Integrate hydrostatic balance for pH with BC of pH(z=0)=0
277           CALL CALC_PH(          CALL CALC_PH(
278       I       bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,       I      bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKp1,
279       U       pH,       U      pH,
280       I       myThid )       I      myThid )
281    
282          ENDDO ! K          ENDDO ! K
283    
284    C--     Initial boundary condition on barotropic divergence integral
285            DO j=1-OLy,sNy+OLy
286             DO i=1-OLx,sNx+OLx
287              cg2d_b(i,j,bi,bj) = 0. _d 0
288             ENDDO
289            ENDDO
290    
291          DO K = Nz, 1, -1          DO K = Nz, 1, -1
292           kM1  =max(1,k-1)   ! Points to level above k (=k-1)           kM1  =max(1,k-1)   ! Points to level above k (=k-1)
293           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above           kUp  =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above
# Line 231  C--      Get temporary terms used by ten Line 303  C--      Get temporary terms used by ten
303       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,       O        xA,yA,uTrans,vTrans,wTrans,maskC,maskUp,
304       I        myThid)       I        myThid)
305    
306  C--      Calculate accelerations in the momentum equations  C--      Calculate the total vertical diffusivity
307           CALL CALC_MOM_RHS(           CALL CALC_DIFFUSIVITY(
308       I        bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I        bi,bj,iMin,iMax,jMin,jMax,K,
309       I        xA,yA,uTrans,vTrans,wTrans,maskC,       I        maskC,maskUp,KapGM,K33,
310       I        pH,       O        KappaZT,
      U        aTerm,xTerm,cTerm,mTerm,pTerm,  
      U        fZon, fMer, fVerU, fVerV,  
311       I        myThid)       I        myThid)
312    
313    
314    C--      Calculate accelerations in the momentum equations
315             IF ( momStepping ) THEN
316              CALL CALC_MOM_RHS(
317         I         bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
318         I         xA,yA,uTrans,vTrans,wTrans,maskC,
319         I         pH,
320         U         aTerm,xTerm,cTerm,mTerm,pTerm,
321         U         fZon, fMer, fVerU, fVerV,
322         I         myThid)
323             ENDIF
324    
325  C--      Calculate active tracer tendencies  C--      Calculate active tracer tendencies
326           CALL CALC_GT(           IF ( tempStepping ) THEN
327       I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,            CALL CALC_GT(
328       I        xA,yA,uTrans,vTrans,wTrans,maskUp,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
329       I        K13,K23,K33,KapGM,       I         xA,yA,uTrans,vTrans,wTrans,maskUp,
330       U        aTerm,xTerm,fZon,fMer,fVerT,       I         K13,K23,KappaZT,KapGM,
331       I        myThid)       U         aTerm,xTerm,fZon,fMer,fVerT,
332         I         myThid)
333             ENDIF
334  Cdbg     CALL CALC_GS(  Cdbg     CALL CALC_GS(
335  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,  Cdbg I        bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
336  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,  Cdbg I        xA,yA,uTrans,vTrans,wTrans,maskUp,
# Line 254  Cdbg I        K13,K23,K33,KapGM, Line 338  Cdbg I        K13,K23,K33,KapGM,
338  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,  Cdbg U        aTerm,xTerm,fZon,fMer,fVerS,
339  Cdbg I        myThid)  Cdbg I        myThid)
340    
341    C--      Prediction step (step forward all model variables)
342             CALL TIMESTEP(
343         I       bi,bj,iMin,iMax,jMin,jMax,K,
344         I       myThid)
345    
346    C--      Diagnose barotropic divergence of predicted fields
347             CALL DIV_G(
348         I       bi,bj,iMin,iMax,jMin,jMax,K,
349         I       xA,yA,
350         I       myThid)
351    
352          ENDDO ! K          ENDDO ! K
353    
354    C--     Implicit diffusion
355            IF (implicitDiffusion) THEN
356             CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
357         I                  KappaZT,
358         I                  myThid )
359            ENDIF
360    
361         ENDDO         ENDDO
362        ENDDO        ENDDO
363    
364  !dbg  write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)        write(0,*) 'dynamics: pS',minval(cg2d_x),maxval(cg2d_x)
365  !dbg  write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: U',minval(uVel(1:sNx,1:sNy,:,:,:)),
366  !dbg &                         maxval(uVel(1:sNx,1:sNy,:,:,:))       &                         maxval(uVel(1:sNx,1:sNy,:,:,:))
367  !dbg  write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: V',minval(vVel(1:sNx,1:sNy,:,:,:)),
368  !dbg &                         maxval(vVel(1:sNx,1:sNy,:,:,:))       &                         maxval(vVel(1:sNx,1:sNy,:,:,:))
369  !dbg  write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
370  !dbg &                         maxval(gT(1:sNx,1:sNy,:,:,:))       &                         maxval(K13(1:sNx,1:sNy,:))
371  !dbg  write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
372  !dbg &                         maxval(Theta(1:sNx,1:sNy,:,:,:))       &                         maxval(K23(1:sNx,1:sNy,:))
373  !dbg  write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),        write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
374  !dbg &                          maxval(pH/(Gravity*Rhonil))       &                         maxval(K33(1:sNx,1:sNy,:))
375          write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),
376         &                         maxval(gT(1:sNx,1:sNy,:,:,:))
377          write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),
378         &                         maxval(Theta(1:sNx,1:sNy,:,:,:))
379          write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),
380         &                          maxval(pH/(Gravity*Rhonil))
381    
382        RETURN        RETURN
383        END        END

Legend:
Removed from v.1.6  
changed lines
  Added in v.1.13

  ViewVC Help
Powered by ViewVC 1.1.22