/[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.11 by adcroft, Mon Jun 1 20:36:13 1998 UTC revision 1.12 by adcroft, Mon Jun 1 22:27:14 1998 UTC
# Line 91  C                          into fVerTerm Line 91  C                          into fVerTerm
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
# Line 106  C       pressure such that the resulting Line 108  C       pressure such that the resulting
108  C       with the free-surface evolution or the rigid-lid:  C       with the free-surface evolution or the rigid-lid:
109  C         U[n] = U* + dt x d/dx P  C         U[n] = U* + dt x d/dx P
110  C         V[n] = V* + dt x d/dy P  C         V[n] = V* + dt x d/dy P
 C       With implicit diffusion, the tracers must also be "finalized"  
 C         (1 + dt * K * d_zz) theta[n] = theta*  
 C         (1 + dt * K * d_zz) salt[n] = salt*  
111  C  C
112  C       "Calculation of Gs"  C       "Calculation of Gs"
113  C       ===================  C       ===================
# Line 122  C         Gv[n] = Gv( u[n], v[n], w, rho Line 121  C         Gv[n] = Gv( u[n], v[n], w, rho
121  C         Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... )  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, ... )  C         Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... )
123  C  C
124  C       "Time-stepping" or "Predicition"  C       "Time-stepping" or "Prediction"
125  C       ================================  C       ================================
126  C       The models variables are stepped forward with the appropriate  C       The models variables are stepped forward with the appropriate
127  C       time-stepping scheme (currently we use Adams-Bashforth II)  C       time-stepping scheme (currently we use Adams-Bashforth II)
# Line 137  C         U* = U[n] + dt x ( 3/2 Gu[n] - Line 136  C         U* = U[n] + dt x ( 3/2 Gu[n] -
136  C         V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )  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] )  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] )  C         salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
139  C       or with implicit diffusion  C       With implicit diffusion:
140  C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )  C         theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
 C  
141  C         salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )  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---  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 167  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
# Line 303  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 the total vertical diffusivity
307             CALL CALC_DIFFUSIVITY(
308         I        bi,bj,iMin,iMax,jMin,jMax,K,
309         I        maskC,maskUp,KapGM,K33,
310         O        KappaZT,
311         I        myThid)
312    
313    
314  C--      Calculate accelerations in the momentum equations  C--      Calculate accelerations in the momentum equations
315           IF ( momStepping ) THEN           IF ( momStepping ) THEN
316            CALL CALC_MOM_RHS(            CALL CALC_MOM_RHS(
# Line 319  C--      Calculate active tracer tendenc Line 327  C--      Calculate active tracer tendenc
327            CALL CALC_GT(            CALL CALC_GT(
328       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,       I         bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown,
329       I         xA,yA,uTrans,vTrans,wTrans,maskUp,       I         xA,yA,uTrans,vTrans,wTrans,maskUp,
330       I         K13,K23,K33,KapGM,       I         K13,K23,KappaZT,KapGM,
331       U         aTerm,xTerm,fZon,fMer,fVerT,       U         aTerm,xTerm,fZon,fMer,fVerT,
332       I         myThid)       I         myThid)
333           ENDIF           ENDIF
# Line 342  C--      Diagnose barotropic divergence Line 350  C--      Diagnose barotropic divergence
350       I       myThid)       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: K13',minval(K13(1:sNx,1:sNy,:)),        write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)),
370  !dbg &                         maxval(K13(1:sNx,1:sNy,:))       &                         maxval(K13(1:sNx,1:sNy,:))
371  !dbg  write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),        write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)),
372  !dbg &                         maxval(K23(1:sNx,1:sNy,:))       &                         maxval(K23(1:sNx,1:sNy,:))
373  !dbg  write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),        write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)),
374  !dbg &                         maxval(K33(1:sNx,1:sNy,:))       &                         maxval(K33(1:sNx,1:sNy,:))
375  !dbg  write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: gT',minval(gT(1:sNx,1:sNy,:,:,:)),
376  !dbg &                         maxval(gT(1:sNx,1:sNy,:,:,:))       &                         maxval(gT(1:sNx,1:sNy,:,:,:))
377  !dbg  write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),        write(0,*) 'dynamics: T',minval(Theta(1:sNx,1:sNy,:,:,:)),
378  !dbg &                         maxval(Theta(1:sNx,1:sNy,:,:,:))       &                         maxval(Theta(1:sNx,1:sNy,:,:,:))
379  !dbg  write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),        write(0,*) 'dynamics: pH',minval(pH/(Gravity*Rhonil)),
380  !dbg &                          maxval(pH/(Gravity*Rhonil))       &                          maxval(pH/(Gravity*Rhonil))
381    
382        RETURN        RETURN
383        END        END

Legend:
Removed from v.1.11  
changed lines
  Added in v.1.12

  ViewVC Help
Powered by ViewVC 1.1.22