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 |
_RL Large_SlopeSqr |
59 |
PARAMETER(Small_Number=1.D-20) |
60 |
PARAMETER(Small_Taper=1.D+03) |
61 |
PARAMETER(Large_SlopeSqr=1.D+08) |
62 |
|
63 |
_RL gradSmod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
64 |
_RL dSigmaDrLtd(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
65 |
_RL dRdSigmaLtd(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
66 |
_RL tmpfld(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
67 |
_RL f1,Smod,f2,Rnondim,Cspd,Lrho |
68 |
_RL maxSlopeSqr |
69 |
_RL fpi |
70 |
PARAMETER(fpi=3.141592653589793047592d0) |
71 |
INTEGER i,j |
72 |
|
73 |
#ifdef ALLOW_AUTODIFF_TAMC |
74 |
act1 = bi - myBxLo(myThid) |
75 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
76 |
act2 = bj - myByLo(myThid) |
77 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
78 |
act3 = myThid - 1 |
79 |
max3 = nTx*nTy |
80 |
act4 = ikey_dynamics - 1 |
81 |
ikey = (act1 + 1) + act2*max1 |
82 |
& + act3*max1*max2 |
83 |
& + act4*max1*max2*max3 |
84 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
85 |
|
86 |
DO j=1-Oly+1,sNy+Oly-1 |
87 |
DO i=1-Olx+1,sNx+Olx-1 |
88 |
gradSmod(i,j) = 0. _d 0 |
89 |
dSigmaDrLtd(i,j) = 0. _d 0 |
90 |
tmpfld(i,j) = 0. _d 0 |
91 |
ENDDO |
92 |
ENDDO |
93 |
|
94 |
IF (GM_taper_scheme.EQ.'orig' .OR. |
95 |
& GM_taper_scheme.EQ.'clipping') THEN |
96 |
|
97 |
#if (defined (GM_TAPER_ORIG_CLIPPING)) |
98 |
|
99 |
C- Original implementation in mitgcmuv |
100 |
C (this turns out to be the same as Cox slope clipping) |
101 |
|
102 |
C- Cox 1987 "Slope clipping" |
103 |
DO j=1-Oly+1,sNy+Oly-1 |
104 |
DO i=1-Olx+1,sNx+Olx-1 |
105 |
tmpfld(i,j) = DSigmaDX(i,j)*DSigmaDX(i,j) + |
106 |
& DSigmaDY(i,j)*DSigmaDY(i,j) |
107 |
IF ( tmpfld(i,j) .EQ. 0. ) THEN |
108 |
gradSmod(i,j) = 0. _d 0 |
109 |
ELSE |
110 |
gradSmod(i,j) = sqrt( tmpfld(i,j) ) |
111 |
ENDIF |
112 |
ENDDO |
113 |
ENDDO |
114 |
|
115 |
#ifdef ALLOW_AUTODIFF_TAMC |
116 |
cnostore CADJ STORE gradSmod(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
117 |
cnostore CADJ STORE dSigmaDrReal(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
118 |
#endif |
119 |
|
120 |
DO j=1-Oly+1,sNy+Oly-1 |
121 |
DO i=1-Olx+1,sNx+Olx-1 |
122 |
IF (gradSmod(i,j) .NE. 0.) THEN |
123 |
dSigmaDrLtd(i,j) = -gradSmod(i,j)*GM_rMaxSlope |
124 |
IF ( dSigmaDrReal(i,j) .GE. dSigmaDrLtd(i,j) ) |
125 |
& dSigmaDrReal(i,j) = dSigmaDrLtd(i,j) |
126 |
ENDIF |
127 |
ENDDO |
128 |
ENDDO |
129 |
|
130 |
#ifdef ALLOW_AUTODIFF_TAMC |
131 |
cnostore CADJ STORE slopeX(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
132 |
cnostore CADJ STORE slopeY(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
133 |
cnostore CADJ STORE dSigmaDrReal(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
134 |
#endif |
135 |
|
136 |
DO j=1-Oly+1,sNy+Oly-1 |
137 |
DO i=1-Olx+1,sNx+Olx-1 |
138 |
IF (gradSmod(i,j) .EQ. 0.) THEN |
139 |
SlopeX(i,j) = 0. _d 0 |
140 |
SlopeY(i,j) = 0. _d 0 |
141 |
ELSE |
142 |
dRdSigmaLtd(i,j) = 1./( dSigmaDrReal(i,j) ) |
143 |
SlopeX(i,j)=-DSigmaDX(i,j)*dRdSigmaLtd(i,j) |
144 |
SlopeY(i,j)=-DSigmaDY(i,j)*dRdSigmaLtd(i,j) |
145 |
ENDIF |
146 |
ENDDO |
147 |
ENDDO |
148 |
|
149 |
#ifdef ALLOW_AUTODIFF_TAMC |
150 |
cnostore CADJ STORE slopeX(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
151 |
cnostore CADJ STORE slopeY(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
152 |
#endif |
153 |
|
154 |
DO j=1-Oly+1,sNy+Oly-1 |
155 |
DO i=1-Olx+1,sNx+Olx-1 |
156 |
SlopeSqr(i,j)=SlopeX(i,j)*SlopeX(i,j) |
157 |
& +SlopeY(i,j)*SlopeY(i,j) |
158 |
taperFct(i,j)=1. _d 0 |
159 |
ENDDO |
160 |
ENDDO |
161 |
|
162 |
#else |
163 |
|
164 |
STOP 'Need to set runtime flag GM_taper_scheme correctly' |
165 |
|
166 |
#endif /* GM_TAPER_ORIG_CLIPPING */ |
167 |
|
168 |
ELSE IF (GM_taper_scheme.EQ.'ac02') THEN |
169 |
|
170 |
#if (defined (GM_TAPER_AC02)) |
171 |
|
172 |
maxSlopeSqr = GM_maxSlope*GM_maxSlope |
173 |
DO j=1-Oly+1,sNy+Oly-1 |
174 |
DO i=1-Olx+1,sNx+Olx-1 |
175 |
dRdSigmaLtd(i,j)= |
176 |
& dSigmaDx(i,j)*dSigmaDx(i,j) |
177 |
& + dSigmaDy(i,j)*dSigmaDy(i,j) |
178 |
& + dSigmaDrReal(i,j)*dSigmaDrReal(i,j) |
179 |
taperFct(i,j) = 1. _d 0 |
180 |
c |
181 |
IF (dRdSigmaLtd(i,j).NE.0.) then |
182 |
dRdSigmaLtd(i,j)=1. _d 0 |
183 |
& / ( dRdSigmaLtd(i,j) ) |
184 |
SlopeSqr(i,j)=(dSigmaDx(i,j)*dSigmaDx(i,j) |
185 |
& +dSigmaDy(i,j)*dSigmaDy(i,j))*dRdSigmaLtd(i,j) |
186 |
SlopeX(i,j)=-dSigmaDx(i,j) |
187 |
& *dRdSigmaLtd(i,j)*dSigmaDrReal(i,j) |
188 |
SlopeY(i,j)=-dSigmaDy(i,j) |
189 |
& *dRdSigmaLtd(i,j)*dSigmaDrReal(i,j) |
190 |
cph T11(i,j)=(dSigmaDrReal(i,j)**2)*dRdSigmaLtd(i,j) |
191 |
ENDIF |
192 |
cph-- this part doesn't adjoint well |
193 |
cph IF ( SlopeSqr(i,j) .GT. maxSlopeSqr .AND. |
194 |
cph & SlopeSqr(i,j) .LT. Large_SlopeSqr ) THEN |
195 |
cph taperFct(i,j) = maxSlopeSqr/SlopeSqr(i,j) |
196 |
cph ELSE IF ( SlopeSqr(i,j) .GT. Large_SlopeSqr ) THEN |
197 |
cph taperFct(i,j) = 0. _d 0 |
198 |
cph ENDIF |
199 |
ENDDO |
200 |
ENDDO |
201 |
|
202 |
#else |
203 |
|
204 |
STOP 'Need to set runtime flag GM_taper_scheme correctly' |
205 |
|
206 |
#endif /* GM_TAPER_AC02 */ |
207 |
|
208 |
ELSE |
209 |
|
210 |
#if (defined (GM_TAPER_REST)) |
211 |
|
212 |
C---------------------------------------------------------------------- |
213 |
|
214 |
C- Compute the slope, no clipping, but avoid reverse slope in negatively |
215 |
C stratified (Sigma_Z > 0) region : |
216 |
|
217 |
#ifdef ALLOW_AUTODIFF_TAMC |
218 |
cnostore CADJ STORE dSigmaDrReal(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
219 |
#endif |
220 |
|
221 |
DO j=1-Oly+1,sNy+Oly-1 |
222 |
DO i=1-Olx+1,sNx+Olx-1 |
223 |
IF ( dSigmaDrReal(i,j) .NE. 0. ) THEN |
224 |
IF (dSigmaDrReal(i,j).GE.(-Small_Number)) |
225 |
& dSigmaDrReal(i,j) = -Small_Number |
226 |
ENDIF |
227 |
ENDDO |
228 |
ENDDO |
229 |
|
230 |
#ifdef ALLOW_AUTODIFF_TAMC |
231 |
cnostore CADJ STORE dSigmaDx(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
232 |
cnostore CADJ STORE dSigmaDy(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
233 |
cnostore CADJ STORE dSigmaDrReal(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
234 |
#endif |
235 |
|
236 |
DO j=1-Oly+1,sNy+Oly-1 |
237 |
DO i=1-Olx+1,sNx+Olx-1 |
238 |
IF ( dSigmaDrReal(i,j) .EQ. 0. ) THEN |
239 |
IF ( dSigmaDx(i,j) .NE. 0. ) THEN |
240 |
SlopeX(i,j) = SIGN(Small_taper,dSigmaDx(i,j)) |
241 |
ELSE |
242 |
SlopeX(i,j) = 0. _d 0 |
243 |
ENDIF |
244 |
IF ( dSigmaDy(i,j) .NE. 0. ) THEN |
245 |
SlopeY(i,j) = SIGN(Small_taper,dSigmaDy(i,j)) |
246 |
ELSE |
247 |
SlopeY(i,j) = 0. _d 0 |
248 |
ENDIF |
249 |
ELSE |
250 |
SlopeX(i,j) = -dSigmaDx(i,j)/dSigmaDrReal(i,j) |
251 |
SlopeY(i,j) = -dSigmaDy(i,j)/dSigmaDrReal(i,j) |
252 |
ENDIF |
253 |
ENDDO |
254 |
ENDDO |
255 |
|
256 |
#ifdef ALLOW_AUTODIFF_TAMC |
257 |
cnostore CADJ STORE slopeX(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
258 |
cnostore CADJ STORE slopeY(:,:) = comlev1_bibj, key=ikey, byte=isbyte |
259 |
#endif |
260 |
|
261 |
DO j=1-Oly+1,sNy+Oly-1 |
262 |
DO i=1-Olx+1,sNx+Olx-1 |
263 |
SlopeSqr(i,j) = SlopeX(i,j)*SlopeX(i,j) |
264 |
& +SlopeY(i,j)*SlopeY(i,j) |
265 |
taperFct(i,j) = 1. _d 0 |
266 |
IF ( SlopeSqr(i,j) .GT. Large_SlopeSqr ) THEN |
267 |
slopeSqr(i,j) = Large_SlopeSqr |
268 |
taperFct(i,j) = 0. _d 0 |
269 |
END IF |
270 |
ENDDO |
271 |
ENDDO |
272 |
|
273 |
C- Compute the tapering function for the GM+Redi tensor : |
274 |
|
275 |
IF (GM_taper_scheme.EQ.'linear') THEN |
276 |
|
277 |
C- Simplest adiabatic tapering = Smax/Slope (linear) |
278 |
maxSlopeSqr = GM_maxSlope*GM_maxSlope |
279 |
DO j=1-Oly+1,sNy+Oly-1 |
280 |
DO i=1-Olx+1,sNx+Olx-1 |
281 |
|
282 |
IF ( SlopeSqr(i,j) .EQ. 0. ) THEN |
283 |
taperFct(i,j) = 1. _d 0 |
284 |
ELSE IF ( SlopeSqr(i,j) .GT. maxSlopeSqr .AND. |
285 |
& SlopeSqr(i,j) .LT. Large_SlopeSqr ) 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 .AND. |
302 |
& SlopeSqr(i,j) .LT. Large_SlopeSqr ) THEN |
303 |
taperFct(i,j) = maxSlopeSqr/SlopeSqr(i,j) |
304 |
ENDIF |
305 |
|
306 |
ENDDO |
307 |
ENDDO |
308 |
|
309 |
ELSEIF (GM_taper_scheme.EQ.'dm95') THEN |
310 |
|
311 |
C- Danabasoglu and McWilliams, J. Clim. 1995 |
312 |
DO j=1-Oly+1,sNy+Oly-1 |
313 |
DO i=1-Olx+1,sNx+Olx-1 |
314 |
|
315 |
IF ( SlopeSqr(i,j) .EQ. 0. ) THEN |
316 |
taperFct(i,j) = 1. _d 0 |
317 |
ELSE IF ( SlopeSqr(i,j) .LT. Large_SlopeSqr ) THEN |
318 |
Smod=sqrt(SlopeSqr(i,j)) |
319 |
taperFct(i,j)=op5*( 1. _d 0 + tanh( (GM_Scrit-Smod)/GM_Sd )) |
320 |
ENDIF |
321 |
ENDDO |
322 |
ENDDO |
323 |
|
324 |
ELSEIF (GM_taper_scheme.EQ.'ldd97') THEN |
325 |
|
326 |
C- Large, Danabasoglu and Doney, JPO 1997 |
327 |
DO j=1-Oly+1,sNy+Oly-1 |
328 |
DO i=1-Olx+1,sNx+Olx-1 |
329 |
|
330 |
IF (SlopeSqr(i,j) .EQ. 0.) THEN |
331 |
taperFct(i,j) = 1. _d 0 |
332 |
ELSE IF ( SlopeSqr(i,j) .LT. Large_SlopeSqr ) THEN |
333 |
Smod=sqrt(SlopeSqr(i,j)) |
334 |
f1=op5*( 1. _d 0 + tanh( (GM_Scrit-Smod)/GM_Sd )) |
335 |
Cspd=2. _d 3 |
336 |
Lrho=100. _d 3 |
337 |
if (FCori(i,j,bi,bj).NE.0.) Lrho=Cspd/abs(Fcori(i,j,bi,bj)) |
338 |
Lrho=min(Lrho , 100. _d 3) |
339 |
Lrho=max(Lrho , 15. _d 3) |
340 |
Rnondim=depthZ/(Lrho*Smod) |
341 |
f2=op5*( 1. _d 0 + sin( fpi*(Rnondim-op5))) |
342 |
taperFct(i,j)=f1*f2 |
343 |
ENDIF |
344 |
|
345 |
ENDDO |
346 |
ENDDO |
347 |
|
348 |
ELSEIF (GM_taper_scheme.NE.' ') THEN |
349 |
STOP 'GMREDI_SLOPE_LIMIT: Bad GM_taper_scheme' |
350 |
ENDIF |
351 |
|
352 |
#endif |
353 |
|
354 |
ENDIF |
355 |
|
356 |
|
357 |
#endif /* ALLOW_GMREDI */ |
358 |
|
359 |
RETURN |
360 |
END |