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

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

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

revision 1.24 by cnh, Sun Feb 4 14:38:45 2001 UTC revision 1.48 by jmc, Sun Sep 23 17:13:13 2007 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2  C $Name$  C $Name$
3    
4    #include "PACKAGES_CONFIG.h"
5  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
6    
7  CStartOfInterFace  CBOP
8        SUBROUTINE CALC_GS(  C     !ROUTINE: CALC_GS
9    C     !INTERFACE:
10          SUBROUTINE CALC_GS(
11       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12       I           xA,yA,uTrans,vTrans,rTrans,maskup,maskC,       I           xA, yA, maskUp, uFld, vFld, wFld,
13         I           uTrans, vTrans, rTrans, rTransKp1,
14       I           KappaRS,       I           KappaRS,
15       U           fVerS,       U           fVerS,
16       I           myCurrentTime, myThid )       I           myTime,myIter,myThid )
17  C     /==========================================================\  C     !DESCRIPTION: \bv
18  C     | SUBROUTINE CALC_GS                                       |  C     *==========================================================*
19  C     | o Calculate the salt tendency terms.                     |  C     | SUBROUTINE CALC_GS
20  C     |==========================================================|  C     | o Calculate the salt tendency terms.
21  C     | A procedure called EXTERNAL_FORCING_S is called from     |  C     *==========================================================*
22  C     | here. These procedures can be used to add per problem    |  C     | A procedure called EXTERNAL_FORCING_S is called from
23  C     | E-P  flux source terms.                                  |  C     | here. These procedures can be used to add per problem
24  C     | Note: Although it is slightly counter-intuitive the      |  C     | E-P  flux source terms.
25  C     |       EXTERNAL_FORCING routine is not the place to put   |  C     | Note: Although it is slightly counter-intuitive the
26  C     |       file I/O. Instead files that are required to       |  C     |       EXTERNAL_FORCING routine is not the place to put
27  C     |       calculate the external source terms are generally  |  C     |       file I/O. Instead files that are required to
28  C     |       read during the model main loop. This makes the    |  C     |       calculate the external source terms are generally
29  C     |       logisitics of multi-processing simpler and also    |  C     |       read during the model main loop. This makes the
30  C     |       makes the adjoint generation simpler. It also      |  C     |       logisitics of multi-processing simpler and also
31  C     |       allows for I/O to overlap computation where that   |  C     |       makes the adjoint generation simpler. It also
32  C     |       is supported by hardware.                          |  C     |       allows for I/O to overlap computation where that
33  C     | Aside from the problem specific term the code here       |  C     |       is supported by hardware.
34  C     | forms the tendency terms due to advection and mixing     |  C     | Aside from the problem specific term the code here
35  C     | The baseline implementation here uses a centered         |  C     | forms the tendency terms due to advection and mixing
36  C     | difference form for the advection term and a tensorial   |  C     | The baseline implementation here uses a centered
37  C     | divergence of a flux form for the diffusive term. The    |  C     | difference form for the advection term and a tensorial
38  C     | diffusive term is formulated so that isopycnal mixing and|  C     | divergence of a flux form for the diffusive term. The
39  C     | GM-style subgrid-scale terms can be incorporated b simply|  C     | diffusive term is formulated so that isopycnal mixing and
40  C     | setting the diffusion tensor terms appropriately.        |  C     | GM-style subgrid-scale terms can be incorporated b simply
41  C     \==========================================================/  C     | setting the diffusion tensor terms appropriately.
42        IMPLICIT NONE  C     *==========================================================*
43    C     \ev
44    
45    C     !USES:
46          IMPLICIT NONE
47  C     == GLobal variables ==  C     == GLobal variables ==
48  #include "SIZE.h"  #include "SIZE.h"
49  #include "DYNVARS.h"  #include "DYNVARS.h"
50  #include "EEPARAMS.h"  #include "EEPARAMS.h"
51  #include "PARAMS.h"  #include "PARAMS.h"
52  #include "GRID.h"  #ifdef ALLOW_GENERIC_ADVDIFF
53  #include "FFIELDS.h"  #include "GAD.h"
54    #endif
55    #ifdef ALLOW_AUTODIFF_TAMC
56    # include "tamc.h"
57    # include "tamc_keys.h"
58    #endif
59    
60    C     !INPUT/OUTPUT PARAMETERS:
61  C     == Routine arguments ==  C     == Routine arguments ==
62  C     fVerS   - Flux of salt (S) in the vertical  C     bi, bj,   :: tile indices
63  C               direction at the upper(U) and lower(D) faces of a cell.  C     iMin,iMax, jMin,jMax :: Range of points for which calculation
64  C     maskUp  - Land mask used to denote base of the domain.  C                                       results will be set.
65  C     maskC   - Land mask for salt cells (used in TOP_LAYER only)  C     k         :: vertical index
66  C     xA      - Tracer cell face area normal to X  C     kM1       :: =k-1 for k>1, =1 for k=1
67  C     yA      - Tracer cell face area normal to X  C     kUp       :: index into 2 1/2D array, toggles between 1|2
68  C     uTrans  - Zonal volume transport through cell face  C     kDown     :: index into 2 1/2D array, toggles between 2|1
69  C     vTrans  - Meridional volume transport through cell face  C     xA        :: Tracer cell face area normal to X
70  C     rTrans  - Vertical volume transport through cell face  C     yA        :: Tracer cell face area normal to X
71  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation  C     maskUp    :: Land mask used to denote base of the domain.
72  C                                      results will be set.  C     uFld,vFld :: Local copy of horizontal velocity field
73  C     myThid - Instance number for this innvocation of CALC_GT  C     wFld      :: Local copy of vertical velocity field
74        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  C     uTrans    :: Zonal volume transport through cell face
75        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     vTrans    :: Meridional volume transport through cell face
76        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     rTrans    ::   Vertical volume transport at interface k
77        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     rTransKp1 :: Vertical volume transport at inteface k+1
78        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     KappaRS   :: Vertical diffusion for Salinity
79        _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     fVerS     :: Flux of salt (S) in the vertical direction
80        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C                  at the upper(U) and lower(D) faces of a cell.
81        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     myTime    :: current time
82        _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  C     myIter    :: current iteration number
83        INTEGER k,kUp,kDown,kM1  C     myThid    :: my Thread Id. number
84    
85        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
86        _RL     myCurrentTime        INTEGER k,kUp,kDown,kM1
87          _RS xA     (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88          _RS yA     (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89          _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90          _RL uFld   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91          _RL vFld   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92          _RL wFld   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93          _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94          _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95          _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96          _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97          _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98          _RL fVerS  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
99          _RL     myTime
100          INTEGER myIter
101        INTEGER myThid        INTEGER myThid
102  CEndOfInterface  CEOP
103    
104    #ifdef ALLOW_GENERIC_ADVDIFF
105    C     === Local variables ===
106          LOGICAL calcAdvection
107          INTEGER iterNb
108    #ifdef ALLOW_ADAMSBASHFORTH_3
109          INTEGER m1, m2
110    #endif
111    
112  C     == Local variables ==  #ifdef ALLOW_AUTODIFF_TAMC
113  C     I, J, K - Loop counters            act1 = bi - myBxLo(myThid)
114        INTEGER i,j            max1 = myBxHi(myThid) - myBxLo(myThid) + 1
115        LOGICAL TOP_LAYER            act2 = bj - myByLo(myThid)
116        _RL afFacS, dfFacS            max2 = myByHi(myThid) - myByLo(myThid) + 1
117        _RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)            act3 = myThid - 1
118        _RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)            max3 = nTx*nTy
119        _RL df4   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)            act4 = ikey_dynamics - 1
120        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)            itdkey = (act1 + 1) + act2*max1
121        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)       &                      + act3*max1*max2
122        _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)       &                      + act4*max1*max2*max3
123        _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)            kkey = (itdkey-1)*Nr + k
124    #endif /* ALLOW_AUTODIFF_TAMC */
125    
126  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF_TAMC
127  C--   only the kUp part of fverS is set in this subroutine  C--   only the kUp part of fverS is set in this subroutine
128  C--   the kDown is still required  C--   the kDown is still required
   
129        fVerS(1,1,kDown) = fVerS(1,1,kDown)        fVerS(1,1,kDown) = fVerS(1,1,kDown)
130        DO j=1-OLy,sNy+OLy  # ifdef NONLIN_FRSURF
131         DO i=1-OLx,sNx+OLx  CADJ STORE fVerS(:,:,:) =
132          fZon(i,j)      = 0.0  CADJ &     comlev1_bibj_k, key=kkey, byte=isbyte
133          fMer(i,j)      = 0.0  CADJ STORE gsNm1(:,:,k,bi,bj) =
134          fVerS(i,j,kUp) = 0.0  CADJ &     comlev1_bibj_k, key=kkey, byte=isbyte
135         ENDDO  # endif
136        ENDDO  #endif
137  #endif  
138    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
139        afFacS = 1. _d 0  
140        dfFacS = 1. _d 0        calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
141        TOP_LAYER = K .EQ. 1        iterNb = myIter
142          IF (staggerTimeStep) iterNb = myIter - 1
143  C---  Calculate advective and diffusive fluxes between cells.  
144    #ifdef ALLOW_ADAMSBASHFORTH_3
145  #ifdef INCLUDE_T_DIFFUSION_CODE          m1 = 1 + MOD(iterNb+1,2)
146  C     o Zonal tracer gradient          m2 = 1 + MOD( iterNb ,2)
147        DO j=1-Oly,sNy+Oly          CALL GAD_CALC_RHS(
148         DO i=1-Olx+1,sNx+Olx       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
149          dSdx(i,j) = _recip_dxC(i,j,bi,bj)*       I           xA, yA, maskUp, uFld, vFld, wFld,
150       &  (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))       I           uTrans, vTrans, rTrans, rTransKp1,
151         ENDDO       I           diffKhS, diffK4S, KappaRS,
152        ENDDO       I           gsNm(1-Olx,1-Oly,1,1,1,m2), salt,
153  C     o Meridional tracer gradient       I           GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
154        DO j=1-Oly+1,sNy+Oly       I           calcAdvection, saltImplVertAdv, AdamsBashforth_S,
155         DO i=1-Olx,sNx+Olx       I           useGMRedi, useKPP,
156          dSdy(i,j) = _recip_dyC(i,j,bi,bj)*       U           fVerS, gS,
157       &  (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))       I           myTime, myIter, myThid )
158         ENDDO  #else /* ALLOW_ADAMSBASHFORTH_3 */
159        ENDDO          CALL GAD_CALC_RHS(
160         I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
161  C--   del^2 of S, needed for bi-harmonic (del^4) term       I           xA, yA, maskUp, uFld, vFld, wFld,
162        IF (diffK4S .NE. 0.) THEN       I           uTrans, vTrans, rTrans, rTransKp1,
163         DO j=1-Oly+1,sNy+Oly-1       I           diffKhS, diffK4S, KappaRS, gsNm1, salt,
164          DO i=1-Olx+1,sNx+Olx-1       I           GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
165           df4(i,j)= _recip_hFacC(i,j,k,bi,bj)       I           calcAdvection, saltImplVertAdv, AdamsBashforth_S,
166       &             *recip_drF(k)/_rA(i,j,bi,bj)       I           useGMRedi, useKPP,
167       &            *(       U           fVerS, gS,
168       &             +( xA(i+1,j)*dSdx(i+1,j)-xA(i,j)*dSdx(i,j) )       I           myTime, myIter, myThid )
169       &             +( yA(i,j+1)*dSdy(i,j+1)-yA(i,j)*dSdy(i,j) )  #endif /* ALLOW_ADAMSBASHFORTH_3 */
170       &             )  
171          ENDDO  C--   External salinity forcing term(s) inside Adams-Bashforth:
172         ENDDO        IF ( saltForcing .AND. tracForcingOutAB.NE.1 )
173        ENDIF       & CALL EXTERNAL_FORCING_S(
174  #endif       I     iMin,iMax,jMin,jMax,bi,bj,k,
175         I     myTime,myThid)
176    
177  C--   Zonal flux (fZon is at west face of "salt" cell)        IF ( AdamsBashforthGs ) THEN
178  C     Advective component of zonal flux  #ifdef ALLOW_ADAMSBASHFORTH_3
179        DO j=jMin,jMax          CALL ADAMS_BASHFORTH3(
180         DO i=iMin,iMax       I                        bi, bj, k,
181          af(i,j) =       U                        gS, gsNm,
182       &   uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0       I                        saltStartAB, iterNb, myThid )
183         ENDDO  #else
184        ENDDO          CALL ADAMS_BASHFORTH2(
185  C     o Diffusive component of zonal flux       I                        bi, bj, k,
186        DO j=jMin,jMax       U                        gS, gsNm1,
187         DO i=iMin,iMax       I                        iterNb, myThid )
         df(i,j) = -diffKhS*xA(i,j)*dSdx(i,j)  
        ENDDO  
       ENDDO  
 #ifdef ALLOW_GMREDI  
       IF (useGMRedi) CALL GMREDI_XTRANSPORT(  
      I     iMin,iMax,jMin,jMax,bi,bj,K,  
      I     xA,salt,  
      U     df,  
      I     myThid)  
 #endif  
 C     o Add the bi-harmonic contribution  
       IF (diffK4S .NE. 0.) THEN  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          df(i,j) = df(i,j) + xA(i,j)*  
      &    diffK4S*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)  
         ENDDO  
        ENDDO  
       ENDIF  
 C     Net zonal flux  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)  
        ENDDO  
       ENDDO  
   
 C--   Meridional flux (fMer is at south face of "salt" cell)  
 C     Advective component of meridional flux  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
 C       Advective component of meridional flux  
         af(i,j) =  
      &   vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0  
        ENDDO  
       ENDDO  
 C     Diffusive component of meridional flux  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         df(i,j) = -diffKhS*yA(i,j)*dSdy(i,j)  
        ENDDO  
       ENDDO  
 #ifdef ALLOW_GMREDI  
       IF (useGMRedi) CALL GMREDI_YTRANSPORT(  
      I     iMin,iMax,jMin,jMax,bi,bj,K,  
      I     yA,salt,  
      U     df,  
      I     myThid)  
188  #endif  #endif
 C     o Add the bi-harmonic contribution  
       IF (diffK4S .NE. 0.) THEN  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          df(i,j) = df(i,j) + yA(i,j)*  
      &    diffK4S*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)  
         ENDDO  
        ENDDO  
189        ENDIF        ENDIF
190    
191  C     Net meridional flux  C--   External salinity forcing term(s) outside Adams-Bashforth:
192        DO j=jMin,jMax        IF ( saltForcing .AND. tracForcingOutAB.EQ.1 )
193         DO i=iMin,iMax       & CALL EXTERNAL_FORCING_S(
194          fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)       I     iMin,iMax,jMin,jMax,bi,bj,k,
195         ENDDO       I     myTime,myThid)
       ENDDO  
   
 C--   Vertical flux (fVerS) above  
 C     Advective component of vertical flux  
 C     Note: For K=1 then KM1=1 this gives a barZ(T) = T  
 C     (this plays the role of the free-surface correction)  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         af(i,j) =  
      &   rTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j,kM1,bi,bj))*0.5 _d 0  
        ENDDO  
       ENDDO  
 C     o Diffusive component of vertical flux  
 C     Note: For K=1 then KM1=1 and this gives a dS/dr = 0 upper  
 C           boundary condition.  
       IF (implicitDiffusion) THEN  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          df(i,j) = 0.  
         ENDDO  
        ENDDO  
       ELSE  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          df(i,j) = - _rA(i,j,bi,bj)*(  
      &    KappaRS(i,j,k)*recip_drC(k)  
      &    *(salt(i,j,kM1,bi,bj)-salt(i,j,k,bi,bj))*rkFac  
      &    )  
         ENDDO  
        ENDDO  
       ENDIF  
   
 #ifdef ALLOW_GMREDI  
       IF (useGMRedi) CALL GMREDI_RTRANSPORT(  
      I     iMin,iMax,jMin,jMax,bi,bj,K,  
      I     maskUp,salt,  
      U     df,  
      I     myThid)  
 #endif  
196    
197  #ifdef ALLOW_KPP  #ifdef NONLIN_FRSURF
198  C--   Add non-local KPP transport term (ghat) to diffusive salt flux.        IF (nonlinFreeSurf.GT.0) THEN
199        IF (useKPP) CALL KPP_TRANSPORT_S(          CALL FREESURF_RESCALE_G(
200       I     iMin,iMax,jMin,jMax,bi,bj,k,km1,       I                          bi, bj, k,
201       I     maskC,KappaRS,       U                          gS,
202       U     df )       I                          myThid )
203            IF ( AdamsBashforthGs ) THEN
204    #ifdef ALLOW_ADAMSBASHFORTH_3
205            CALL FREESURF_RESCALE_G(
206         I                          bi, bj, k,
207         U                          gsNm(1-OLx,1-OLy,1,1,1,1),
208         I                          myThid )
209            CALL FREESURF_RESCALE_G(
210         I                          bi, bj, k,
211         U                          gsNm(1-OLx,1-OLy,1,1,1,2),
212         I                          myThid )
213    #else
214            CALL FREESURF_RESCALE_G(
215         I                          bi, bj, k,
216         U                          gsNm1,
217         I                          myThid )
218  #endif  #endif
219            ENDIF
 C     Net vertical flux  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         fVerS(i,j,kUp) = ( afFacS*af(i,j)+  dfFacS*df(i,j) )*maskUp(i,j)  
        ENDDO  
       ENDDO  
       IF ( TOP_LAYER ) THEN  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          fVerS(i,j,kUp) = afFacS*af(i,j)*freeSurfFac  
         ENDDO  
        ENDDO  
220        ENDIF        ENDIF
221    #endif /* NONLIN_FRSURF */
222    
223  C--   Tendency is minus divergence of the fluxes.  #endif /* ALLOW_GENERIC_ADVDIFF */
 C     Note. Tendency terms will only be correct for range  
 C           i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points  
 C           will contain valid floating point numbers but  
 C           they are not algorithmically correct. These points  
 C           are not used.  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
 #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)  
 #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)  
         gS(i,j,k,bi,bj)=  
      &   -_recip_VolS1(i,j,k,bi,bj)  
      &    _recip_VolS2(i,j,k,bi,bj)  
      &   *(  
      &    +( fZon(i+1,j)-fZon(i,j) )  
      &    +( fMer(i,j+1)-fMer(i,j) )  
      &    +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac  
      &    )  
        ENDDO  
       ENDDO  
   
 C--   External forcing term(s)  
       CALL EXTERNAL_FORCING_S(  
      I     iMin,iMax,jMin,jMax,bi,bj,k,  
      I     maskC,  
      I     myCurrentTime,myThid)  
224    
225        RETURN        RETURN
226        END        END

Legend:
Removed from v.1.24  
changed lines
  Added in v.1.48

  ViewVC Help
Powered by ViewVC 1.1.22