/[MITgcm]/MITgcm/pkg/gmredi/gmredi_slope_psi.F
ViewVC logotype

Diff of /MITgcm/pkg/gmredi/gmredi_slope_psi.F

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

revision 1.2 by heimbach, Fri Feb 15 21:28:07 2002 UTC revision 1.8 by jmc, Fri Nov 4 01:31:36 2005 UTC
# Line 4  C $Name$ Line 4  C $Name$
4  #include "GMREDI_OPTIONS.h"  #include "GMREDI_OPTIONS.h"
5    
6  CStartOfInterface  CStartOfInterface
7        SUBROUTINE GMREDI_SLOPE_PSI_B(        SUBROUTINE GMREDI_SLOPE_PSI(
      I             dSigmaDrW,dSigmaDrS,  
      I             depthZ,  
      U             SlopeX, SlopeY,  
8       O             taperX, taperY,       O             taperX, taperY,
9         U             SlopeX, SlopeY,
10         U             dSigmaDrW,dSigmaDrS,
11         I             LrhoW, LrhoS, depthZ, K,
12       I             bi,bj, myThid )       I             bi,bj, myThid )
13  C     /==========================================================\  C     /==========================================================\
14  C     | SUBROUTINE GMREDI_SLOPE_PSI_B                            |  C     | SUBROUTINE GMREDI_SLOPE_PSI                              |
15  C     | o Calculate slopes for use in GM/Redi tensor             |  C     | o Calculate slopes for use in GM/Redi tensor             |
16  C     |==========================================================|  C     |==========================================================|
17  C     | On entry:                                                |  C     | On entry:                                                |
18  C     |            dSigmaDrW conatins the d/dz Sigma          |  C     |            dSigmaDrW,S  contains the d/dz Sigma          |
19  C     |            SlopeX/Y     contains X/Y gradients of sigma  |  C     |            SlopeX/Y     contains X/Y gradients of sigma  |
20  C     |            depthZ       conatins the height (m) of level |  C     |            depthZ       contains the height (m) of level |
21  C     | On exit:                                                 |  C     | On exit:                                                 |
22  C     |            dSigmaDrW conatins the effective dSig/dz   |  C     |            dSigmaDrW,S  contains the effective dSig/dz   |
23  C     |            SlopeX/Y     contains X/Y slopes              |  C     |            SlopeX/Y     contains X/Y slopes              |
24  C     |            taperFct     contains tapering funct. value ; |  C     |            taperFct     contains tapering funct. value ; |
25  C     |                         = 1 when using no tapering       |  C     |                         = 1 when using no tapering       |
# Line 32  C     == Global variables == Line 32  C     == Global variables ==
32  #include "GMREDI.h"  #include "GMREDI.h"
33  #include "PARAMS.h"  #include "PARAMS.h"
34    
35    #ifdef ALLOW_AUTODIFF_TAMC
36    #include "tamc.h"
37    #include "tamc_keys.h"
38    #endif /* ALLOW_AUTODIFF_TAMC */
39    
40  C     == Routine arguments ==  C     == Routine arguments ==
41  C  C
42          _RL taperX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
43          _RL taperY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
44        _RL SlopeX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL SlopeX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
45        _RL SlopeY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL SlopeY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
46        _RL dSigmaDrW(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL dSigmaDrW(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
47        _RL dSigmaDrS(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL dSigmaDrS(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
48        _RL taperX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL LrhoW(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
49        _RL taperY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL LrhoS(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
50        _RL depthZ        _RL depthZ
51        INTEGER bi,bj,myThid        INTEGER K,bi,bj,myThid
52  CEndOfInterface  CEndOfInterface
53    
54  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
55    #ifdef GM_BOLUS_ADVEC
56    
57  C     == Local variables ==  C     == Local variables ==
58        _RL Small_Number        _RL dSigmaDrLtd(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
59        PARAMETER(Small_Number=1.D-12)        _RL f1,Smod,f2,Rnondim
       _RL gradSmod,f1,Smod,f2,Rnondim,Cspd,Lrho  
       _RL dSigmaDrLtd, dRdSigmaLtd  
60        _RL maxSlopeSqr        _RL maxSlopeSqr
61          _RL slopeCutoff
62        _RL fpi        _RL fpi
63        PARAMETER(fpi=3.141592653589793047592d0)        PARAMETER(fpi=3.141592653589793047592d0)
64        INTEGER i,j        INTEGER i,j
65    
66          slopeCutoff = SQRT( GM_slopeSqCutoff )
67    
68    #ifdef ALLOW_AUTODIFF_TAMC
69          act1 = bi - myBxLo(myThid)
70          max1 = myBxHi(myThid) - myBxLo(myThid) + 1
71          act2 = bj - myByLo(myThid)
72          max2 = myByHi(myThid) - myByLo(myThid) + 1
73          act3 = myThid - 1
74          max3 = nTx*nTy
75          act4 = ikey_dynamics - 1
76          igmkey = (act1 + 1) + act2*max1
77         &                    + act3*max1*max2
78         &                    + act4*max1*max2*max3
79          kkey = (igmkey-1)*Nr + k
80    #endif /* ALLOW_AUTODIFF_TAMC */
81    
82        IF (GM_taper_scheme.EQ.'orig' .OR.        IF (GM_taper_scheme.EQ.'orig' .OR.
83       &    GM_taper_scheme.EQ.'clipping') THEN       &    GM_taper_scheme.EQ.'clipping') THEN
84    
85    #ifdef GM_EXCLUDE_CLIPPING
86    
87            STOP 'Need to compile without "#define GM_EXCLUDE_CLIPPING"'
88    
89    #else  /* GM_EXCLUDE_CLIPPING */
90    
91  C-      Original implementation in mitgcmuv  C-      Original implementation in mitgcmuv
92  C       (this turns out to be the same as Cox slope clipping)  C       (this turns out to be the same as Cox slope clipping)
93    
94  C-      Cox 1987 "Slope clipping"  C-- X-comp
         DO j=1-Oly+1,sNy+Oly-1  
          DO i=1-Olx+1,sNx+Olx-1  
95    
96  c         gradSmod=SlopeX(i,j)*SlopeX(i,j)  #ifdef ALLOW_AUTODIFF_TAMC
97  c    &            +SlopeY(i,j)*SlopeY(i,j)          DO j=1-Oly,sNy+Oly
98  c         if (gradSmod .NE. 0.) gradSmod=sqrt(gradSmod)           DO i=1-Olx+1,sNx+Olx
99            gradSmod=abs(SlopeX(i,j))            dSigmaDrLtd(i,j) = 0. _d 0
100             ENDDO
101            dSigmaDrLtd = -(Small_Number+gradSmod*GM_rMaxSlope)          ENDDO
102            IF (dSigmaDrW(i,j).GE.dSigmaDrLtd)  #endif /* ALLOW_AUTODIFF_TAMC */
103       &        dSigmaDrW(i,j) = dSigmaDrLtd  
104    C-      Cox 1987 "Slope clipping"
105            DO j=1-Oly,sNy+Oly
106             DO i=1-Olx+1,sNx+Olx
107              dSigmaDrLtd(i,j) = -(GM_Small_Number+
108         &     ABS(SlopeX(i,j))*GM_rMaxSlope)
109             ENDDO
110            ENDDO
111    #ifdef ALLOW_AUTODIFF_TAMC
112    CADJ STORE dSigmaDrLtd(:,:)  = comlev1_bibj_k, key=kkey, byte=isbyte
113    CADJ STORE dSigmaDrW(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
114    #endif
115            DO j=1-Oly,sNy+Oly
116             DO i=1-Olx+1,sNx+Olx
117              IF (dSigmaDrW(i,j).GE.dSigmaDrLtd(i,j))
118         &        dSigmaDrW(i,j) = dSigmaDrLtd(i,j)
119             ENDDO
120            ENDDO
121    #ifdef ALLOW_AUTODIFF_TAMC
122    CADJ STORE dSigmaDrW(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
123    #endif
124            DO j=1-Oly,sNy+Oly
125             DO i=1-Olx+1,sNx+Olx
126            SlopeX(i,j) = -SlopeX(i,j)/dSigmaDrW(i,j)            SlopeX(i,j) = -SlopeX(i,j)/dSigmaDrW(i,j)
127              taperX(i,j)=1. _d 0
128             ENDDO
129            ENDDO
130    
131            gradSmod=abs(SlopeY(i,j))  C-- Y-comp
           dSigmaDrLtd = -(Small_Number+gradSmod*GM_rMaxSlope)  
           IF (dSigmaDrS(i,j).GE.dSigmaDrLtd)  
      &        dSigmaDrS(i,j) = dSigmaDrLtd  
           SlopeY(i,j) = -SlopeY(i,j)/dSigmaDrS(i,j)  
132    
133            taperX(i,j)=1. _d 0  #ifdef ALLOW_AUTODIFF_TAMC
134            DO j=1-Oly+1,sNy+Oly
135             DO i=1-Olx,sNx+Olx
136              dSigmaDrLtd(i,j) = 0. _d 0
137             ENDDO
138            ENDDO
139    #endif /* ALLOW_AUTODIFF_TAMC */
140            DO j=1-Oly+1,sNy+Oly
141             DO i=1-Olx,sNx+Olx
142              dSigmaDrLtd(i,j) = -(GM_Small_Number+
143         &     ABS(SlopeY(i,j))*GM_rMaxSlope)
144             ENDDO
145            ENDDO
146    #ifdef ALLOW_AUTODIFF_TAMC
147    CADJ STORE dSigmaDrLtd(:,:)  = comlev1_bibj_k, key=kkey, byte=isbyte
148    CADJ STORE dSigmaDrS(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
149    #endif
150            DO j=1-Oly+1,sNy+Oly
151             DO i=1-Olx,sNx+Olx
152              IF (dSigmaDrS(i,j).GE.dSigmaDrLtd(i,j))
153         &        dSigmaDrS(i,j) = dSigmaDrLtd(i,j)
154             ENDDO
155            ENDDO
156    #ifdef ALLOW_AUTODIFF_TAMC
157    CADJ STORE dSigmaDrS(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
158    #endif
159            DO j=1-Oly+1,sNy+Oly
160             DO i=1-Olx,sNx+Olx
161              SlopeY(i,j) = -SlopeY(i,j)/dSigmaDrS(i,j)
162            taperY(i,j)=1. _d 0            taperY(i,j)=1. _d 0
   
163           ENDDO           ENDDO
164          ENDDO          ENDDO
165    
166    #endif /* GM_EXCLUDE_CLIPPING */
167    
168        ELSE        ELSE
169    
170    #ifdef GM_EXCLUDE_TAPERING
171    
172            STOP 'Need to compile without "#define GM_EXCLUDE_TAPERING"'
173    
174    #else  /* GM_EXCLUDE_TAPERING */
175    
176    #ifdef ALLOW_AUTODIFF_TAMC
177    CADJ STORE slopeX(:,:)       = comlev1_bibj_k, key=kkey, byte=isbyte
178    CADJ STORE dSigmaDrW(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
179    #endif
180    
181  C- Compute the slope, no clipping, but avoid reverse slope in negatively  C- Compute the slope, no clipping, but avoid reverse slope in negatively
182  C                                  stratified (Sigma_Z > 0) region :  C                                  stratified (Sigma_Z > 0) region :
183          DO j=1-Oly+1,sNy+Oly-1          DO j=1-Oly,sNy+Oly
184           DO i=1-Olx+1,sNx+Olx-1           DO i=1-Olx+1,sNx+Olx
185              IF (dSigmaDrW(i,j).GE.-GM_Small_Number)
186            dSigmaDrLtd = -Small_Number       &        dSigmaDrW(i,j) = -GM_Small_Number
187            IF (dSigmaDrW(i,j).GE.dSigmaDrLtd)           ENDDO
188       &        dSigmaDrW(i,j) = dSigmaDrLtd          ENDDO
189            dRdSigmaLtd = 1./dSigmaDrW(i,j)  #ifdef ALLOW_AUTODIFF_TAMC
190    CADJ STORE dsigmadrW(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
191    #endif
192            DO j=1-Oly,sNy+Oly
193             DO i=1-Olx+1,sNx+Olx
194            SlopeX(i,j) = -SlopeX(i,j)/dSigmaDrW(i,j)            SlopeX(i,j) = -SlopeX(i,j)/dSigmaDrW(i,j)
195              taperX(i,j)= 1. _d 0
196             ENDDO
197            ENDDO
198    #ifdef ALLOW_AUTODIFF_TAMC
199    CADJ STORE slopex(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
200    #endif
201            DO j=1-Oly,sNy+Oly
202             DO i=1-Olx+1,sNx+Olx
203              IF ( ABS(SlopeX(i,j)) .GE. slopeCutoff ) THEN
204                 SlopeX(i,j) = SIGN(slopeCutoff,SlopeX(i,j))
205                 taperX(i,j) = 0. _d 0
206              ENDIF
207             ENDDO
208            ENDDO
209    
210            dSigmaDrLtd = -Small_Number  #ifdef ALLOW_AUTODIFF_TAMC
211            IF (dSigmaDrS(i,j).GE.dSigmaDrLtd)  CADJ STORE slopeY(:,:)       = comlev1_bibj_k, key=kkey, byte=isbyte
212       &        dSigmaDrS(i,j) = dSigmaDrLtd  CADJ STORE dSigmaDrS(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
213    #endif
214    
215            DO j=1-Oly+1,sNy+Oly
216             DO i=1-Olx,sNx+Olx
217              IF (dSigmaDrS(i,j).GE.-GM_Small_Number)
218         &        dSigmaDrS(i,j) = -GM_Small_Number
219             ENDDO
220            ENDDO
221    #ifdef ALLOW_AUTODIFF_TAMC
222    CADJ STORE dsigmadrS(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
223    #endif
224            DO j=1-Oly+1,sNy+Oly
225             DO i=1-Olx,sNx+Olx
226            SlopeY(i,j) = -SlopeY(i,j)/dSigmaDrS(i,j)            SlopeY(i,j) = -SlopeY(i,j)/dSigmaDrS(i,j)
   
 c         SlopeSqr(i,j)=SlopeX(i,j)*SlopeX(i,j)  
 c    &                 +SlopeY(i,j)*SlopeY(i,j)  
   
           taperX(i,j)=1. _d 0  
227            taperY(i,j)=1. _d 0            taperY(i,j)=1. _d 0
228             ENDDO
229            ENDDO
230    #ifdef ALLOW_AUTODIFF_TAMC
231    CADJ STORE slopey(:,:)    = comlev1_bibj_k, key=kkey, byte=isbyte
232    #endif
233            DO j=1-Oly+1,sNy+Oly
234             DO i=1-Olx,sNx+Olx
235              IF ( ABS(SlopeY(i,j)) .GE. slopeCutoff ) THEN
236                 SlopeY(i,j) = SIGN(slopeCutoff,SlopeY(i,j))
237                 taperY(i,j) = 0. _d 0
238              ENDIF
239           ENDDO           ENDDO
240          ENDDO          ENDDO
241    
242  C- Compute the tapering function for the GM+Redi tensor :  C- Compute the tapering function for the GM+Redi tensor :
243    
244    #ifdef ALLOW_AUTODIFF_TAMC
245    CADJ STORE slopeX(:,:)       = comlev1_bibj_k, key=kkey, byte=isbyte
246    CADJ STORE slopeY(:,:)       = comlev1_bibj_k, key=kkey, byte=isbyte
247    #endif
248    
249         IF (GM_taper_scheme.EQ.'linear') THEN         IF (GM_taper_scheme.EQ.'linear') THEN
250    
251  C-      Simplest adiabatic tapering = Smax/Slope (linear)  C-      Simplest adiabatic tapering = Smax/Slope (linear)
252          DO j=1-Oly+1,sNy+Oly-1          DO j=1-Oly,sNy+Oly
253           DO i=1-Olx+1,sNx+Olx-1           DO i=1-Olx+1,sNx+Olx
254              Smod = ABS(SlopeX(i,j))
255            IF (abs(SlopeX(i,j)).GT.GM_maxSlope)            IF ( Smod .GT. GM_maxSlope .AND.
256       &      taperX(i,j)=GM_maxSlope/abs(SlopeX(i,j))       &           Smod .LT. slopeCutoff )
257            IF (abs(SlopeY(i,j)).GT.GM_maxSlope)       &           taperX(i,j)=GM_maxSlope/(Smod+GM_Small_Number)
258       &      taperY(i,j)=GM_maxSlope/abs(SlopeY(i,j))           ENDDO
259            ENDDO
260            DO j=1-Oly+1,sNy+Oly
261             DO i=1-Olx,sNx+Olx
262              Smod = ABS(SlopeY(i,j))
263              IF ( Smod .GT. GM_maxSlope .AND.
264         &           Smod .LT. slopeCutoff )
265         &           taperY(i,j)=GM_maxSlope/(Smod+GM_Small_Number)
266           ENDDO           ENDDO
267          ENDDO          ENDDO
268    
# Line 135  C-      Simplest adiabatic tapering = Sm Line 270  C-      Simplest adiabatic tapering = Sm
270    
271  C-      Gerdes, Koberle and Willebrand, Clim. Dyn. 1991  C-      Gerdes, Koberle and Willebrand, Clim. Dyn. 1991
272          maxSlopeSqr = GM_maxSlope*GM_maxSlope          maxSlopeSqr = GM_maxSlope*GM_maxSlope
273          DO j=1-Oly+1,sNy+Oly-1          DO j=1-Oly,sNy+Oly
274           DO i=1-Olx+1,sNx+Olx-1           DO i=1-Olx+1,sNx+Olx
275              IF ( ABS(SlopeX(i,j)) .GT. GM_maxSlope .AND.
276            IF (abs(SlopeX(i,j)).GT.GM_maxSlope)       &           ABS(SlopeX(i,j)) .LT. slopeCutoff )
277       &      taperX(i,j)=maxSlopeSqr/(SlopeX(i,j)*SlopeX(i,j))       &           taperX(i,j)=maxSlopeSqr/
278            IF (abs(SlopeY(i,j)).GT.GM_maxSlope)       &           ( SlopeX(i,j)*SlopeX(i,j) + GM_Small_Number )
279       &      taperY(i,j)=maxSlopeSqr/(SlopeY(i,j)*SlopeY(i,j))           ENDDO
280            ENDDO
281            DO j=1-Oly+1,sNy+Oly
282             DO i=1-Olx,sNx+Olx
283              IF ( ABS(SlopeY(i,j)) .GT. GM_maxSlope .AND.
284         &           ABS(SlopeY(i,j)) .LT. slopeCutoff )
285         &           taperY(i,j)=maxSlopeSqr/
286         &           ( SlopeY(i,j)*SlopeY(i,j) + GM_Small_Number )
287           ENDDO           ENDDO
288          ENDDO          ENDDO
289    
290         ELSEIF (GM_taper_scheme.EQ.'dm95') THEN         ELSEIF (GM_taper_scheme.EQ.'dm95') THEN
291    
292  C-      Danabasoglu and McWilliams, J. Clim. 1995  C-      Danabasoglu and McWilliams, J. Clim. 1995
293          DO j=1-Oly+1,sNy+Oly-1          DO j=1-Oly,sNy+Oly
294           DO i=1-Olx+1,sNx+Olx-1           DO i=1-Olx+1,sNx+Olx
295              Smod = ABS(SlopeX(i,j))
296            Smod = abs(SlopeX(i,j))            taperX(i,j)=op5*( 1. _d 0 + TANH( (GM_Scrit-Smod)/GM_Sd ))
297            taperX(i,j)=0.5*(1.+tanh( (GM_Scrit-Smod)/GM_Sd ))           ENDDO
298            Smod = abs(SlopeY(i,j))          ENDDO
299            taperY(i,j)=0.5*(1.+tanh( (GM_Scrit-Smod)/GM_Sd ))          DO j=1-Oly+1,sNy+Oly
300             DO i=1-Olx,sNx+Olx
301              Smod = ABS(SlopeY(i,j))
302              taperY(i,j)=op5*( 1. _d 0 + TANH( (GM_Scrit-Smod)/GM_Sd ))
303           ENDDO           ENDDO
304          ENDDO          ENDDO
305    
306         ELSEIF (GM_taper_scheme.EQ.'ldd97') THEN         ELSEIF (GM_taper_scheme.EQ.'ldd97') THEN
307    
308  C-      Large, Danabasoglu and Doney, JPO 1997  C-      Large, Danabasoglu and Doney, JPO 1997
         DO j=1-Oly+1,sNy+Oly-1  
          DO i=1-Olx+1,sNx+Olx-1  
309    
310            Cspd=2.          DO j=1-Oly,sNy+Oly
311            Lrho=100.e3           DO i=1-Olx+1,sNx+Olx
312            if (FCori(i,j,bi,bj).NE.0.) Lrho=Cspd/abs(Fcori(i,j,bi,bj))            Smod = ABS(SlopeX(i,j))
313            Lrho=min(Lrho , 100. _d 3)            IF ( Smod .LT. slopeCutoff ) THEN
314            Lrho=max(Lrho , 15. _d 3)            f1=op5*( 1. _d 0 + TANH( (GM_Scrit-Smod)/GM_Sd ))
315              IF (Smod.NE.0.) THEN
316            Smod = abs(SlopeX(i,j))              Rnondim=depthZ/(LrhoW(i,j)*Smod)
317            f1=0.5*(1.+tanh( (GM_Scrit-Smod)/GM_Sd ))            ELSE
           if (Smod.NE.0.) then  
             Rnondim=depthZ/(Lrho*Smod)  
           else  
318              Rnondim=0.              Rnondim=0.
319            endif            ENDIF
320            f2=0.5*(1.+sin( fpi*(Rnondim-0.5)))            f2=op5*( 1. _d 0 + SIN( fpi*(Rnondim-op5)))
321            taperX(i,j)=f1*f2            taperX(i,j)=f1*f2
322              ENDIF
323             ENDDO
324            ENDDO
325    
326            Smod = abs(SlopeY(i,j))          DO j=1-Oly+1,sNy+Oly
327            f1=0.5*(1.+tanh( (GM_Scrit-Smod)/GM_Sd ))           DO i=1-Olx,sNx+Olx
328            if (Smod.NE.0.) then            Smod = ABS(SlopeY(i,j))
329              Rnondim=depthZ/(Lrho*Smod)            IF ( Smod .LT. slopeCutoff ) THEN
330            else            f1=op5*( 1. _d 0 + TANH( (GM_Scrit-Smod)/GM_Sd ))
331              IF (Smod.NE.0.) THEN
332                Rnondim=depthZ/(LrhoS(i,j)*Smod)
333              ELSE
334              Rnondim=0.              Rnondim=0.
335            endif            ENDIF
336            f2=0.5*(1.+sin( fpi*(Rnondim-0.5)))            f2=op5*( 1. _d 0 + SIN( fpi*(Rnondim-op5)))
337            taperY(i,j)=f1*f2            taperY(i,j)=f1*f2
338              ENDIF
339           ENDDO           ENDDO
340          ENDDO          ENDDO
341    
# Line 199  C-      Large, Danabasoglu and Doney, JP Line 343  C-      Large, Danabasoglu and Doney, JP
343          STOP 'GMREDI_SLOPE_PSI: Bad GM_taper_scheme'          STOP 'GMREDI_SLOPE_PSI: Bad GM_taper_scheme'
344         ENDIF         ENDIF
345    
346        ENDIF  #endif /* GM_EXCLUDE_TAPERING */
347    
348          ENDIF
349    
350    #endif /* BOLUS_ADVEC */
351  #endif /* ALLOW_GMREDI */  #endif /* ALLOW_GMREDI */
352    
353        RETURN        RETURN

Legend:
Removed from v.1.2  
changed lines
  Added in v.1.8

  ViewVC Help
Powered by ViewVC 1.1.22