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

Contents of /MITgcm/pkg/gmredi/gmredi_slope_limit.F

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


Revision 1.14 - (show annotations) (download)
Fri Nov 15 02:57:47 2002 UTC (21 years, 5 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint47e_post, checkpoint47c_post, checkpoint47d_pre, checkpoint47a_post, checkpoint47d_post, branch-exfmods-tag, checkpoint47b_post, checkpoint47f_post, checkpoint47
Branch point for: branch-exfmods-curt
Changes since 1.13: +16 -6 lines
o * "clean" adjoint code (in terms of extensive recomputations)
    can now be obtained for all GMREDI options (i.e. for
    - GM_VISBECK_VARIABLE_K
    - GM_NON_UNITY_DIAGONAL
    - GM_EXTRA_DIAGONAL
    - GM_BOLUS_ADVEC )
  * However, wrong gradient check problem remains unsolved.
  * New CPP options have been introduced for different
    tapering schemes

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

  ViewVC Help
Powered by ViewVC 1.1.22