/[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.15 by cnh, Wed Oct 28 03:11:36 1998 UTC revision 1.47 by jmc, Sun Jun 18 23:22:43 2006 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    C $Name$
3    
4  #include "CPP_EEOPTIONS.h"  #include "PACKAGES_CONFIG.h"
5    #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           K13,K23,KappaRS,KapGM,       I           uTrans, vTrans, rTrans, rTransKp1,
14       U           af,df,fZon,fMer,fVerS,       I           KappaRS,
15       I           myThid )       U           fVerS,
16  C     /==========================================================\       I           myTime,myIter,myThid )
17  C     | SUBROUTINE CALC_GS                                       |  C     !DESCRIPTION: \bv
18  C     | o Calculate the salt tendency terms.                     |  C     *==========================================================*
19  C     |==========================================================|  C     | SUBROUTINE CALC_GS
20  C     | A procedure called EXTERNAL_FORCING_S is called from     |  C     | o Calculate the salt tendency terms.
21  C     | here. These procedures can be used to add per problem    |  C     *==========================================================*
22  C     | E-P  flux source terms.                                  |  C     | A procedure called EXTERNAL_FORCING_S is called from
23  C     | Note: Although it is slightly counter-intuitive the      |  C     | here. These procedures can be used to add per problem
24  C     |       EXTERNAL_FORCING routine is not the place to put   |  C     | E-P  flux source terms.
25  C     |       file I/O. Instead files that are required to       |  C     | Note: Although it is slightly counter-intuitive the
26  C     |       calculate the external source terms are generally  |  C     |       EXTERNAL_FORCING routine is not the place to put
27  C     |       read during the model main loop. This makes the    |  C     |       file I/O. Instead files that are required to
28  C     |       logisitics of multi-processing simpler and also    |  C     |       calculate the external source terms are generally
29  C     |       makes the adjoint generation simpler. It also      |  C     |       read during the model main loop. This makes the
30  C     |       allows for I/O to overlap computation where that   |  C     |       logisitics of multi-processing simpler and also
31  C     |       is supported by hardware.                          |  C     |       makes the adjoint generation simpler. It also
32  C     | Aside from the problem specific term the code here       |  C     |       allows for I/O to overlap computation where that
33  C     | forms the tendency terms due to advection and mixing     |  C     |       is supported by hardware.
34  C     | The baseline implementation here uses a centered         |  C     | Aside from the problem specific term the code here
35  C     | difference form for the advection term and a tensorial   |  C     | forms the tendency terms due to advection and mixing
36  C     | divergence of a flux form for the diffusive term. The    |  C     | The baseline implementation here uses a centered
37  C     | diffusive term is formulated so that isopycnal mixing and|  C     | difference form for the advection term and a tensorial
38  C     | GM-style subgrid-scale terms can be incorporated b simply|  C     | divergence of a flux form for the diffusive term. The
39  C     | setting the diffusion tensor terms appropriately.        |  C     | diffusive term is formulated so that isopycnal mixing and
40  C     \==========================================================/  C     | GM-style subgrid-scale terms can be incorporated b simply
41        IMPLICIT NONE  C     | setting the diffusion tensor terms appropriately.
42    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     fZon    - Work array for flux of temperature in the east-west  C     bi, bj,   :: tile indices
63  C               direction at the west face of a cell.  C     iMin,iMax, jMin,jMax :: Range of points for which calculation
64  C     fMer    - Work array for flux of temperature in the north-south  C                                       results will be set.
65  C               direction at the south face of a cell.  C     k         :: vertical index
66  C     fVerS   - Flux of salt (S) in the vertical  C     kM1       :: =k-1 for k>1, =1 for k=1
67  C               direction at the upper(U) and lower(D) faces of a cell.  C     kUp       :: index into 2 1/2D array, toggles between 1|2
68  C     maskUp  - Land mask used to denote base of the domain.  C     kDown     :: index into 2 1/2D array, toggles between 2|1
69  C     maskC   - Land mask for salt cells (used in TOP_LAYER only)  C     xA        :: Tracer cell face area normal to X
70  C     xA      - Tracer cell face area normal to X  C     yA        :: Tracer cell face area normal to X
71  C     yA      - Tracer cell face area normal to X  C     maskUp    :: Land mask used to denote base of the domain.
72  C     uTrans  - Zonal volume transport through cell face  C     uFld,vFld :: Local copy of horizontal velocity field
73  C     vTrans  - Meridional volume transport through cell face  C     wFld      :: Local copy of vertical velocity field
74  C     wTrans  - Vertical volume transport through cell face  C     uTrans    :: Zonal volume transport through cell face
75  C     af      - Advective flux component work array  C     vTrans    :: Meridional volume transport through cell face
76  C     df      - Diffusive flux component work array  C     rTrans    ::   Vertical volume transport at interface k
77  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation  C     rTransKp1 :: Vertical volume transport at inteface k+1
78  C                                      results will be set.  C     KappaRS   :: Vertical diffusion for Salinity
79  C     myThid - Instance number for this innvocation of CALC_GT  C     fVerS     :: Flux of salt (S) in the vertical direction
80        _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C                  at the upper(U) and lower(D) faces of a cell.
81        _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     myTime    :: current time
82        _RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)  C     myIter    :: current iteration number
83        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  C     myThid    :: my Thread Id. number
84        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL KappaRS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)  
       _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       INTEGER k,kUp,kDown,kM1  
85        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
86          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  C     == Local variables ==  #ifdef ALLOW_GENERIC_ADVDIFF
105  C     I, J, K - Loop counters  C     === Local variables ===
106        INTEGER i,j        LOGICAL calcAdvection
107        LOGICAL TOP_LAYER        INTEGER iterNb
108        _RL afFacS, dfFacS  #ifdef ALLOW_ADAMSBASHFORTH_3
109        _RL dSdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        INTEGER m1, m2
110        _RL dSdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)  #endif
111    
112        afFacS = 1. _d 0  #ifdef ALLOW_AUTODIFF_TAMC
113        dfFacS = 1. _d 0            act1 = bi - myBxLo(myThid)
114        TOP_LAYER = K .EQ. 1            max1 = myBxHi(myThid) - myBxLo(myThid) + 1
115              act2 = bj - myByLo(myThid)
116  C---  Calculate advective and diffusive fluxes between cells.            max2 = myByHi(myThid) - myByLo(myThid) + 1
117              act3 = myThid - 1
118  C--   Zonal flux (fZon is at west face of "salt" cell)            max3 = nTx*nTy
119  C     Advective component of zonal flux            act4 = ikey_dynamics - 1
120        DO j=jMin,jMax            itdkey = (act1 + 1) + act2*max1
121         DO i=iMin,iMax       &                      + act3*max1*max2
122          af(i,j) =       &                      + act4*max1*max2*max3
123       &   uTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i-1,j,k,bi,bj))*0.5 _d 0            kkey = (itdkey-1)*Nr + k
124         ENDDO  #endif /* ALLOW_AUTODIFF_TAMC */
125        ENDDO  
126  C     Zonal tracer gradient  #ifdef ALLOW_AUTODIFF_TAMC
127        DO j=jMin,jMax  C--   only the kUp part of fverS is set in this subroutine
128         DO i=iMin,iMax  C--   the kDown is still required
129          dSdx(i,j) = _recip_dxC(i,j,bi,bj)*        fVerS(1,1,kDown) = fVerS(1,1,kDown)
130       &  (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj))  # ifdef NONLIN_FRSURF
131         ENDDO  CADJ STORE fVerS(:,:,:) =
132        ENDDO  CADJ &     comlev1_bibj_k, key=kkey, byte=isbyte
133  C     Diffusive component of zonal flux  CADJ STORE gsNm1(:,:,k,bi,bj) =
134        DO j=jMin,jMax  CADJ &     comlev1_bibj_k, key=kkey, byte=isbyte
135         DO i=iMin,iMax  # endif
136          df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i-1,j)))*  #endif
137       &            xA(i,j)*dSdx(i,j)  
138         ENDDO  C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
139        ENDDO  
140  C     Net zonal flux        calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
141        DO j=jMin,jMax        iterNb = myIter
142         DO i=iMin,iMax        IF (staggerTimeStep) iterNb = myIter - 1
143          fZon(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)  
144         ENDDO  #ifdef ALLOW_ADAMSBASHFORTH_3
145        ENDDO          m1 = 1 + MOD(iterNb+1,2)
146            m2 = 1 + MOD( iterNb ,2)
147  C--   Meridional flux (fMer is at south face of "salt" cell)          CALL GAD_CALC_RHS(
148  C     Advective component of meridional flux       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
149        DO j=jMin,jMax       I           xA, yA, maskUp, uFld, vFld, wFld,
150         DO i=iMin,iMax       I           uTrans, vTrans, rTrans, rTransKp1,
151  C       Advective component of meridional flux       I           diffKhS, diffK4S, KappaRS,
152          af(i,j) =       I           gsNm(1-Olx,1-Oly,1,1,1,m2), salt,
153       &   vTrans(i,j)*(salt(i,j,k,bi,bj)+salt(i,j-1,k,bi,bj))*0.5 _d 0       I           GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
154         ENDDO       I           calcAdvection, saltImplVertAdv, AdamsBashforth_S,
155        ENDDO       U           fVerS, gS,
156  C     Zonal tracer gradient       I           myTime, myIter, myThid )
157        DO j=jMin,jMax  #else /* ALLOW_ADAMSBASHFORTH_3 */
158         DO i=iMin,iMax          CALL GAD_CALC_RHS(
159          dSdy(i,j) = _recip_dyC(i,j,bi,bj)*       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
160       &  (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))       I           xA, yA, maskUp, uFld, vFld, wFld,
161         ENDDO       I           uTrans, vTrans, rTrans, rTransKp1,
162        ENDDO       I           diffKhS, diffK4S, KappaRS, gsNm1, salt,
163  C     Diffusive component of meridional flux       I           GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
164        DO j=jMin,jMax       I           calcAdvection, saltImplVertAdv, AdamsBashforth_S,
165         DO i=iMin,iMax       U           fVerS, gS,
166          df(i,j) = -(diffKhS+0.5*(KapGM(i,j)+KapGM(i,j-1)))*       I           myTime, myIter, myThid )
167       &            yA(i,j)*dSdy(i,j)  #endif /* ALLOW_ADAMSBASHFORTH_3 */
168         ENDDO  
169        ENDDO  C--   External salinity forcing term(s) inside Adams-Bashforth:
170  C     Net meridional flux        IF ( saltForcing .AND. tracForcingOutAB.NE.1 )
171        DO j=jMin,jMax       & CALL EXTERNAL_FORCING_S(
172         DO i=iMin,iMax       I     iMin,iMax,jMin,jMax,bi,bj,k,
173          fMer(i,j) = afFacS*af(i,j) + dfFacS*df(i,j)       I     myTime,myThid)
174         ENDDO  
175        ENDDO        IF ( AdamsBashforthGs ) THEN
176    #ifdef ALLOW_ADAMSBASHFORTH_3
177  C--   Interpolate terms for Redi/GM scheme          CALL ADAMS_BASHFORTH3(
178        DO j=jMin,jMax       I                        bi, bj, k,
179         DO i=iMin,iMax       U                        gS, gsNm,
180          dSdx(i,j) = 0.5*(       I                        saltStartAB, iterNb, myThid )
181       &   +0.5*(_maskW(i+1,j,k,bi,bj)  #else
182       &         *_recip_dxC(i+1,j,bi,bj)*          CALL ADAMS_BASHFORTH2(
183       &           (salt(i+1,j,k,bi,bj)-salt(i,j,k,bi,bj))       I                        bi, bj, k,
184       &        +_maskW(i,j,k,bi,bj)       U                        gS, gsNm1,
185       &         *_recip_dxC(i,j,bi,bj)*       I                        iterNb, myThid )
186       &           (salt(i,j,k,bi,bj)-salt(i-1,j,k,bi,bj)))  #endif
      &   +0.5*(_maskW(i+1,j,km1,bi,bj)  
      &         *_recip_dxC(i+1,j,bi,bj)*  
      &           (salt(i+1,j,km1,bi,bj)-salt(i,j,km1,bi,bj))  
      &        +_maskW(i,j,km1,bi,bj)  
      &         *_recip_dxC(i,j,bi,bj)*  
      &           (salt(i,j,km1,bi,bj)-salt(i-1,j,km1,bi,bj)))  
      &       )  
        ENDDO  
       ENDDO  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         dSdy(i,j) = 0.5*(  
      &   +0.5*(_maskS(i,j,k,bi,bj)  
      &         *_recip_dyC(i,j,bi,bj)*  
      &           (salt(i,j,k,bi,bj)-salt(i,j-1,k,bi,bj))  
      &        +_maskS(i,j+1,k,bi,bj)  
      &         *_recip_dyC(i,j+1,bi,bj)*  
      &           (salt(i,j+1,k,bi,bj)-salt(i,j,k,bi,bj)))  
      &   +0.5*(_maskS(i,j,km1,bi,bj)  
      &         *_recip_dyC(i,j,bi,bj)*  
      &           (salt(i,j,km1,bi,bj)-salt(i,j-1,km1,bi,bj))  
      &        +_maskS(i,j+1,km1,bi,bj)  
      &         *_recip_dyC(i,j+1,bi,bj)*  
      &           (salt(i,j+1,km1,bi,bj)-salt(i,j,km1,bi,bj)))  
      &       )  
        ENDDO  
       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     Diffusive component of vertical flux  
 C     Note: For K=1 then KM1=1 this gives a dS/dz = 0 upper  
 C           boundary condition.  
       DO j=jMin,jMax  
        DO i=iMin,iMax  
         df(i,j) = _rA(i,j,bi,bj)*(  
      &   -KapGM(i,j)*K13(i,j,k)*dSdx(i,j)  
      &   -KapGM(i,j)*K23(i,j,k)*dSdy(i,j)  
      &   )  
        ENDDO  
       ENDDO  
       IF (.NOT.implicitDiffusion) THEN  
        DO j=jMin,jMax  
         DO i=iMin,iMax  
          df(i,j) = 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  
 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  
187        ENDIF        ENDIF
188    
189  C--   Tendency is minus divergence of the fluxes.  C--   External salinity forcing term(s) outside Adams-Bashforth:
190  C     Note. Tendency terms will only be correct for range        IF ( saltForcing .AND. tracForcingOutAB.EQ.1 )
191  C           i=iMin+1:iMax-1, j=jMin+1:jMax-1. Edge points       & CALL EXTERNAL_FORCING_S(
192  C           will contain valid floating point numbers but       I     iMin,iMax,jMin,jMax,bi,bj,k,
193  C           they are not algorithmically correct. These points       I     myTime,myThid)
194  C           are not used.  
195        DO j=jMin,jMax  #ifdef NONLIN_FRSURF
196         DO i=iMin,iMax        IF (nonlinFreeSurf.GT.0) THEN
197  #define _recip_VolS1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)          CALL FREESURF_RESCALE_G(
198  #define _recip_VolS2(i,j,k,bi,bj) /_rA(i,j,bi,bj)       I                          bi, bj, k,
199          gS(i,j,k,bi,bj)=       U                          gS,
200       &   -_recip_VolS1(i,j,k,bi,bj)       I                          myThid )
201       &    _recip_VolS2(i,j,k,bi,bj)          IF ( AdamsBashforthGs ) THEN
202       &   *(  #ifdef ALLOW_ADAMSBASHFORTH_3
203       &    +( fZon(i+1,j)-fZon(i,j) )          CALL FREESURF_RESCALE_G(
204       &    +( fMer(i,j+1)-fMer(i,j) )       I                          bi, bj, k,
205       &    +( fVerS(i,j,kUp)-fVerS(i,j,kDown) )*rkFac       U                          gsNm(1-OLx,1-OLy,1,1,1,1),
206       &    )       I                          myThid )
207         ENDDO          CALL FREESURF_RESCALE_G(
208        ENDDO       I                          bi, bj, k,
209         U                          gsNm(1-OLx,1-OLy,1,1,1,2),
210  C--   External P-E forcing term(s)       I                          myThid )
211  C     o Surface relaxation term  #else
212        IF ( TOP_LAYER ) THEN          CALL FREESURF_RESCALE_G(
213         DO j=jMin,jMax       I                          bi, bj, k,
214          DO i=iMin,iMax       U                          gsNm1,
215           gS(i,j,k,bi,bj)=gS(i,j,k,bi,bj)       I                          myThid )
216       &   +maskC(i,j)*(  #endif
217       &   -lambdaSaltClimRelax*(salt(i,j,k,bi,bj)-SSS(i,j,bi,bj))          ENDIF
      &   +EmPmR(i,j,bi,bj) )  
         ENDDO  
        ENDDO  
218        ENDIF        ENDIF
219    #endif /* NONLIN_FRSURF */
220    
221    #endif /* ALLOW_GENERIC_ADVDIFF */
222    
223        RETURN        RETURN
224        END        END

Legend:
Removed from v.1.15  
changed lines
  Added in v.1.47

  ViewVC Help
Powered by ViewVC 1.1.22