/[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.19 - (show annotations) (download)
Tue Jan 21 19:34:13 2003 UTC (21 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint50c_post, checkpoint48e_post, checkpoint50c_pre, checkpoint48i_post, checkpoint51, checkpoint50, checkpoint50d_post, checkpoint50b_pre, checkpoint48b_post, checkpoint51d_post, checkpoint48c_pre, checkpoint48d_pre, checkpoint48d_post, checkpoint48f_post, checkpoint48h_post, checkpoint51b_pre, checkpoint48a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint47j_post, branchpoint-genmake2, checkpoint48c_post, checkpoint51b_post, checkpoint51c_post, checkpoint50g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint50e_post, checkpoint50d_pre, checkpoint51e_post, checkpoint48, checkpoint49, checkpoint51f_pre, checkpoint48g_post, checkpoint50b_post, checkpoint51a_post
Branch point for: branch-genmake2
Changes since 1.18: +9 -7 lines
Yet more changes:
o adgmredi_calc_tensor
  avoiding all recomputation of gmredi_slope_limit
o adgmredi_x/y/rtransport
  added flag for excessive storing to avoid recomp. of
  u/v/rtans, dTdx/y/z
  -> this is not really necessary and very memory-consuming
o adgmredi_slope_psi:
  consistency with gmredi_slope_limit in treatment of GM_slopeSqCutoff
o gmredi_slope_limit
  re-activated full calculation of taperfct for case 'ac02'

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

  ViewVC Help
Powered by ViewVC 1.1.22