/[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.10 by adcroft, Thu May 28 16:19:50 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
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 120  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 130  C     uninitialised but inert locations. Line 177  C     uninitialised but inert locations.
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    
 C--     Boundary condition on hydrostatic pressure is pH(z=0)=0  
         DO j=1-OLy,sNy+OLy  
          DO i=1-OLx,sNx+OLx  
           pH(i,j,1) = 0. _d 0  
           K13(i,j,1) = 0. _d 0  
           K23(i,j,1) = 0. _d 0  
           K33(i,j,1) = 0. _d 0  
           KapGM(i,j) = 0. _d 0  
          ENDDO  
         ENDDO  
   
180  C--     Set up work arrays that need valid initial values  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
# Line 153  C--     Set up work arrays that need val Line 189  C--     Set up work arrays that need val
189            fVerU(i,j,2) = 0. _d 0            fVerU(i,j,2) = 0. _d 0
190            fVerV(i,j,1) = 0. _d 0            fVerV(i,j,1) = 0. _d 0
191            fVerV(i,j,2) = 0. _d 0            fVerV(i,j,2) = 0. _d 0
192              pH(i,j,1) = 0. _d 0
193              K13(i,j,1) = 0. _d 0
194              K23(i,j,1) = 0. _d 0
195              K33(i,j,1) = 0. _d 0
196              KapGM(i,j) = 0. _d 0
197           ENDDO           ENDDO
198          ENDDO          ENDDO
199    
# Line 168  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
# Line 189  C--     Integrate hydrostatic balance fo Line 230  C--     Integrate hydrostatic balance fo
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-1 level  C--     Density of K-1 level (above W(K)) reference to K-1 level
236  copt    CALL FIND_RHO(  copt    CALL FIND_RHO(
# Line 238  C--     Integrate hydrostatic balance fo Line 279  C--     Integrate hydrostatic balance fo
279       U      pH,       U      pH,
280       I      myThid )       I      myThid )
281    
282            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          ENDDO
290    
291          DO K = Nz, 1, -1          DO K = Nz, 1, -1
# Line 255  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 271  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 282  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          ENDDO  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
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.10  
changed lines
  Added in v.1.12

  ViewVC Help
Powered by ViewVC 1.1.22