1 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_dynsolver.F,v 1.15 2007/03/08 11:21:34 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "SEAICE_OPTIONS.h" |
5 |
|
6 |
CStartOfInterface |
7 |
SUBROUTINE SEAICE_DYNSOLVER( myTime, myIter, myThid ) |
8 |
C /==========================================================\ |
9 |
C | SUBROUTINE SEAICE_DYNSOLVER | |
10 |
C | o Ice dynamics using LSR solver | |
11 |
C | Zhang and Hibler, JGR, 102, 8691-8702, 1997 | |
12 |
C | or EVP explicit solver by Hunke and Dukowicz, JPO 27, | |
13 |
C | 1849-1867 (1997) | |
14 |
C |==========================================================| |
15 |
C | written by Martin Losch, March 2006 | |
16 |
C \==========================================================/ |
17 |
IMPLICIT NONE |
18 |
|
19 |
C === Global variables === |
20 |
#include "SIZE.h" |
21 |
#include "EEPARAMS.h" |
22 |
#include "PARAMS.h" |
23 |
#include "GRID.h" |
24 |
#include "SURFACE.h" |
25 |
#include "DYNVARS.h" |
26 |
#include "FFIELDS.h" |
27 |
#include "SEAICE.h" |
28 |
#include "SEAICE_PARAMS.h" |
29 |
#include "SEAICE_FFIELDS.h" |
30 |
|
31 |
#ifdef ALLOW_AUTODIFF_TAMC |
32 |
# include "tamc.h" |
33 |
#endif |
34 |
|
35 |
C === Routine arguments === |
36 |
C myTime - Simulation time |
37 |
C myIter - Simulation timestep number |
38 |
C myThid - Thread no. that called this routine. |
39 |
_RL myTime |
40 |
INTEGER myIter |
41 |
INTEGER myThid |
42 |
CEndOfInterface |
43 |
|
44 |
#ifdef SEAICE_CGRID |
45 |
C === Local variables === |
46 |
C i,j,bi,bj - Loop counters |
47 |
|
48 |
INTEGER i, j, bi, bj |
49 |
_RL RHOICE, RHOAIR, SINWIN, COSWIN |
50 |
_RL PSTAR, AAA |
51 |
_RL U1, V1 |
52 |
_RL phiSurf(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
53 |
|
54 |
C-- FIRST SET UP BASIC CONSTANTS |
55 |
RHOICE = SEAICE_rhoIce |
56 |
RHOAIR = SEAICE_rhoAir |
57 |
PSTAR = SEAICE_strength |
58 |
|
59 |
C-- introduce turning angle (default is zero) |
60 |
SINWIN=SIN(SEAICE_airTurnAngle*deg2rad) |
61 |
COSWIN=COS(SEAICE_airTurnAngle*deg2rad) |
62 |
|
63 |
C-- Compute proxy for geostrophic velocity, |
64 |
DO bj=myByLo(myThid),myByHi(myThid) |
65 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
66 |
DO j=1-Oly,sNy+Oly |
67 |
DO i=1-Oly,sNx+Olx |
68 |
CML GWATX(I,J,bi,bj)=uVel(i,j,KGEO(I,J,bi,bj),bi,bj) |
69 |
CML GWATY(I,J,bi,bj)=vVel(i,j,KGEO(I,J,bi,bj),bi,bj) |
70 |
GWATX(I,J,bi,bj)=uVel(i,j,1,bi,bj) |
71 |
GWATY(I,J,bi,bj)=vVel(i,j,1,bi,bj) |
72 |
ENDDO |
73 |
ENDDO |
74 |
ENDDO |
75 |
ENDDO |
76 |
|
77 |
C-- NOW SET UP MASS PER UNIT AREA AND CORIOLIS TERM |
78 |
DO bj=myByLo(myThid),myByHi(myThid) |
79 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
80 |
DO j=1-Oly+1,sNy+Oly |
81 |
DO i=1-Olx+1,sNx+Olx |
82 |
seaiceMassC(I,J,bi,bj)=RHOICE*HEFF(i,j,1,bi,bj) |
83 |
seaiceMassU(I,J,bi,bj)=RHOICE*HALF*( |
84 |
& HEFF(i,j,1,bi,bj) + HEFF(i-1,j ,1,bi,bj) ) |
85 |
seaiceMassV(I,J,bi,bj)=RHOICE*HALF*( |
86 |
& HEFF(i,j,1,bi,bj) + HEFF(i ,j-1,1,bi,bj) ) |
87 |
ENDDO |
88 |
ENDDO |
89 |
ENDDO |
90 |
ENDDO |
91 |
|
92 |
IF ( SEAICE_maskRHS ) THEN |
93 |
C dynamic masking of areas with no ice |
94 |
DO bj=myByLo(myThid),myByHi(myThid) |
95 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
96 |
DO j=1-Oly+1,sNy+Oly |
97 |
DO i=1-Olx+1,sNx+Olx |
98 |
seaiceMaskU(I,J,bi,bj)=AREA(I,J,1,bi,bj)+AREA(I-1,J,1,bi,bj) |
99 |
IF ( seaiceMaskU(I,J,bi,bj) .GT. 0. _d 0 ) THEN |
100 |
seaiceMaskU(I,J,bi,bj) = 1. _d 0 |
101 |
ELSE |
102 |
seaiceMaskU(I,J,bi,bj) = 0. _d 0 |
103 |
ENDIF |
104 |
seaiceMaskV(I,J,bi,bj)=AREA(I,J,1,bi,bj)+AREA(I,J-1,1,bi,bj) |
105 |
IF ( seaiceMaskV(I,J,bi,bj) .GT. 0. _d 0 ) THEN |
106 |
seaiceMaskV(I,J,bi,bj) = 1. _d 0 |
107 |
ELSE |
108 |
seaiceMaskV(I,J,bi,bj) = 0. _d 0 |
109 |
ENDIF |
110 |
ENDDO |
111 |
ENDDO |
112 |
ENDDO |
113 |
ENDDO |
114 |
CALL EXCH_UV_XY_RL( seaiceMaskU, seaiceMaskV, .FALSE., myThid ) |
115 |
ENDIF |
116 |
|
117 |
C-- NOW SET UP FORCING FIELDS |
118 |
|
119 |
#ifdef ALLOW_ATM_WIND |
120 |
C-- Wind stress is computed on center of C-grid cell and interpolated |
121 |
C to U and V points later |
122 |
C locations from wind on tracer locations |
123 |
DO bj=myByLo(myThid),myByHi(myThid) |
124 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
125 |
DO j=1-Oly,sNy+Oly |
126 |
DO i=1-Olx,sNx+Olx |
127 |
U1=UWIND(I,J,bi,bj) |
128 |
V1=VWIND(I,J,bi,bj) |
129 |
AAA=U1**2+V1**2 |
130 |
IF ( AAA .LE. SEAICE_EPS_SQ ) THEN |
131 |
AAA=SEAICE_EPS |
132 |
ELSE |
133 |
AAA=SQRT(AAA) |
134 |
ENDIF |
135 |
C first ocean surface stress |
136 |
DAIRN(I,J,bi,bj)=RHOAIR*OCEAN_drag |
137 |
& *(2.70 _d 0+0.142 _d 0*AAA+0.0764 _d 0*AAA*AAA) |
138 |
WINDX(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
139 |
& (COSWIN*U1-SIGN(SINWIN, _fCori(I,J,bi,bj))*V1) |
140 |
WINDY(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
141 |
& (SIGN(SINWIN, _fCori(I,J,bi,bj))*U1+COSWIN*V1) |
142 |
C now ice surface stress |
143 |
DAIRN(I,J,bi,bj) = RHOAIR*SEAICE_drag*AAA |
144 |
ENDDO |
145 |
ENDDO |
146 |
C now interpolate to U and V points respectively |
147 |
DO j=1-Oly+1,sNy+Oly |
148 |
DO i=1-Olx+1,sNx+Olx |
149 |
FORCEX0(I,J,bi,bj)=0.5 _d 0 * |
150 |
& ( DAIRN(I ,J,bi,bj)*( |
151 |
& COSWIN*uWind(I ,J,bi,bj) |
152 |
& -SIGN(SINWIN, _fCori(I ,J,bi,bj))*vWind(I ,J,bi,bj) ) |
153 |
& + DAIRN(I-1,J,bi,bj)*( |
154 |
& COSWIN*uWind(I-1,J,bi,bj) |
155 |
& -SIGN(SINWIN, _fCori(I-1,J,bi,bj))*vWind(I-1,J,bi,bj) ) |
156 |
& ) |
157 |
C interpolate to V point |
158 |
FORCEY0(I,J,bi,bj)=0.5 _d 0 * |
159 |
& ( DAIRN(I,J ,bi,bj)*( |
160 |
& SIGN(SINWIN, _fCori(I,J ,bi,bj))*uWind(I,J ,bi,bj) |
161 |
& +COSWIN*vWind(I,J ,bi,bj) ) |
162 |
& + DAIRN(I,J-1,bi,bj)*( |
163 |
& SIGN(SINWIN, _fCori(I,J-1,bi,bj))*uWind(I,J-1,bi,bj) |
164 |
& +COSWIN*vWind(I,J-1,bi,bj) ) |
165 |
& ) |
166 |
ENDDO |
167 |
ENDDO |
168 |
ENDDO |
169 |
ENDDO |
170 |
#else |
171 |
C-- Wind stress is available on U and V points, copy it to seaice variables. |
172 |
DO bj=myByLo(myThid),myByHi(myThid) |
173 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
174 |
DO j=1-Oly,sNy+Oly-1 |
175 |
DO i=1-Olx,sNx+Olx-1 |
176 |
U1=.5*(FU(I,J,bi,bj)+FU(I+1,J,bi,bj)) |
177 |
V1=.5*(FV(I,J,bi,bj)+FV(I,J+1,bi,bj)) |
178 |
C first ocean surface stress |
179 |
WINDX(I,J,bi,bj)= |
180 |
& (COSWIN*U1 |
181 |
& -SIGN(SINWIN, _fCori(I,J,bi,bj))*V1) |
182 |
WINDY(I,J,bi,bj)= |
183 |
& (SIGN(SINWIN, _fCori(I,J,bi,bj))*U1+COSWIN*V1) |
184 |
ENDDO |
185 |
ENDDO |
186 |
C now interpolate to U and V points respectively |
187 |
DO j=1-Oly+1,sNy+Oly-1 |
188 |
DO i=1-Olx+1,sNx+Olx-1 |
189 |
C now ice surface stress |
190 |
DAIRN(I,J,bi,bj) = SEAICE_drag/OCEAN_drag |
191 |
FORCEX0(I,J,bi,bj)=DAIRN(I,J,bi,bj)*( |
192 |
& COSWIN*FU(I,J,bi,bj) |
193 |
& -SIGN(SINWIN, _fCori(I ,J,bi,bj))*0.25 |
194 |
& *(FV(I,J, bi,bj)+FV(I-1,J, bi,bj) |
195 |
& +FV(I,J+1,bi,bj)+FV(I-1,J+1,bi,bj)) |
196 |
& ) |
197 |
C interpolate to V point |
198 |
FORCEY0(I,J,bi,bj)=DAIRN(I,J,bi,bj)*( |
199 |
& SIGN(SINWIN, _fCori(I,J ,bi,bj))*0.25 |
200 |
& *(FU(I,J ,bi,bj)+FU(I+1,J, bi,bj) |
201 |
& +FU(I,J-1,bi,bj)+FU(I+1,J-1,bi,bj)) |
202 |
& +COSWIN*FV(I,J,bi,bj) |
203 |
& ) |
204 |
ENDDO |
205 |
ENDDO |
206 |
ENDDO |
207 |
ENDDO |
208 |
#endif /* ALLOW_ATM_WIND */ |
209 |
|
210 |
DO bj=myByLo(myThid),myByHi(myThid) |
211 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
212 |
C-- Compute surface pressure at z==0: |
213 |
C- use actual sea surface height for tilt computations |
214 |
DO j=1-Oly,sNy+Oly |
215 |
DO i=1-Olx,sNx+Olx |
216 |
phiSurf(i,j) = Bo_surf(i,j,bi,bj)*etaN(i,j,bi,bj) |
217 |
ENDDO |
218 |
ENDDO |
219 |
#ifdef ATMOSPHERIC_LOADING |
220 |
C- add atmospheric loading and Sea-Ice loading |
221 |
IF ( useRealFreshWaterFlux ) THEN |
222 |
DO j=1-Oly,sNy+Oly |
223 |
DO i=1-Olx,sNx+Olx |
224 |
phiSurf(i,j) = phiSurf(i,j) |
225 |
& + ( pload(i,j,bi,bj) |
226 |
& +sIceLoad(i,j,bi,bj)*gravity |
227 |
& )*recip_rhoConst |
228 |
ENDDO |
229 |
ENDDO |
230 |
ELSE |
231 |
DO j=1-Oly,sNy+Oly |
232 |
DO i=1-Olx,sNx+Olx |
233 |
phiSurf(i,j) = phiSurf(i,j) |
234 |
& + pload(i,j,bi,bj)*recip_rhoConst |
235 |
ENDDO |
236 |
ENDDO |
237 |
ENDIF |
238 |
#endif /* ATMOSPHERIC_LOADING */ |
239 |
DO j=1-Oly+1,sNy+Oly |
240 |
DO i=1-Olx+1,sNx+Olx |
241 |
C-- NOW ADD IN TILT |
242 |
FORCEX0(I,J,bi,bj)=FORCEX0(I,J,bi,bj) |
243 |
& -seaiceMassU(I,J,bi,bj)*_recip_dxC(I,J,bi,bj) |
244 |
& *( phiSurf(i,j)-phiSurf(i-1,j) ) |
245 |
FORCEY0(I,J,bi,bj)=FORCEY0(I,J,bi,bj) |
246 |
& -seaiceMassV(I,J,bi,bj)* _recip_dyC(I,J,bi,bj) |
247 |
& *( phiSurf(i,j)-phiSurf(i,j-1) ) |
248 |
C-- NOW SET UP ICE PRESSURE AND VISCOSITIES |
249 |
PRESS0(I,J,bi,bj)=PSTAR*HEFF(I,J,1,bi,bj) |
250 |
& *EXP(-20.0 _d 0*(ONE-AREA(I,J,1,bi,bj))) |
251 |
ZMAX(I,J,bi,bj)=(5.0 _d +12/(2.0 _d +04))*PRESS0(I,J,bi,bj) |
252 |
ZMIN(I,J,bi,bj)=4.0 _d +08 |
253 |
PRESS0(I,J,bi,bj)=PRESS0(I,J,bi,bj)*HEFFM(I,J,bi,bj) |
254 |
ENDDO |
255 |
ENDDO |
256 |
ENDDO |
257 |
ENDDO |
258 |
|
259 |
#ifdef SEAICE_ALLOW_DYNAMICS |
260 |
|
261 |
IF ( SEAICEuseDYNAMICS ) THEN |
262 |
|
263 |
#ifdef ALLOW_AUTODIFF_TAMC |
264 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
265 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
266 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
267 |
|
268 |
#ifdef SEAICE_ALLOW_EVP |
269 |
IF ( SEAICEuseEVP ) THEN |
270 |
CALL SEAICE_EVP( myTime, myIter, myThid ) |
271 |
ELSE |
272 |
C do the original VP solver of Zhang and Hibler (1997), ported to |
273 |
C a C-grid |
274 |
#endif /* SEAICE_ALLOW_EVP */ |
275 |
C NOW DO PREDICTOR TIME STEP |
276 |
DO bj=myByLo(myThid),myByHi(myThid) |
277 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
278 |
DO j=1-OLy,sNy+OLy |
279 |
DO i=1-OLx,sNx+OLx |
280 |
uIce(I,J,2,bi,bj)=uIce(I,J,1,bi,bj) |
281 |
vIce(I,J,2,bi,bj)=vIce(I,J,1,bi,bj) |
282 |
uIceC(I,J,bi,bj)=uIce(I,J,1,bi,bj) |
283 |
vIceC(I,J,bi,bj)=vIce(I,J,1,bi,bj) |
284 |
ENDDO |
285 |
ENDDO |
286 |
ENDDO |
287 |
ENDDO |
288 |
|
289 |
#ifdef ALLOW_AUTODIFF_TAMC |
290 |
# ifdef SEAICE_ALLOW_EVP |
291 |
cphCADJ STORE uice = comlev1, key=ikey_dynamics |
292 |
cphCADJ STORE vice = comlev1, key=ikey_dynamics |
293 |
CADJ STORE uicec = comlev1, key=ikey_dynamics |
294 |
CADJ STORE vicec = comlev1, key=ikey_dynamics |
295 |
# endif |
296 |
#endif |
297 |
C NOW LSR SCHEME (ZHANG-J/HIBLER 1997) |
298 |
CALL SEAICE_LSR( 1, myThid ) |
299 |
#ifdef ALLOW_AUTODIFF_TAMC |
300 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
301 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
302 |
cphCADJ STORE uicec = comlev1, key=ikey_dynamics |
303 |
cphCADJ STORE vicec = comlev1, key=ikey_dynamics |
304 |
#endif |
305 |
|
306 |
C NOW DO MODIFIED EULER STEP |
307 |
DO bj=myByLo(myThid),myByHi(myThid) |
308 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
309 |
DO j=1-OLy,sNy+OLy |
310 |
DO i=1-OLx,sNx+OLx |
311 |
uIce(I,J,1,bi,bj)=HALF*(uIce(I,J,1,bi,bj)+uIce(I,J,2,bi,bj)) |
312 |
vIce(I,J,1,bi,bj)=HALF*(vIce(I,J,1,bi,bj)+vIce(I,J,2,bi,bj)) |
313 |
uIceC(I,J,bi,bj)=uIce(I,J,1,bi,bj) |
314 |
vIceC(I,J,bi,bj)=vIce(I,J,1,bi,bj) |
315 |
ENDDO |
316 |
ENDDO |
317 |
ENDDO |
318 |
ENDDO |
319 |
|
320 |
C NOW LSR SCHEME (ZHANG-J/HIBLER 1997) |
321 |
CALL SEAICE_LSR( 2, myThid ) |
322 |
#ifdef SEAICE_ALLOW_EVP |
323 |
ENDIF |
324 |
#endif /* SEAICE_ALLOW_EVP */ |
325 |
|
326 |
ENDIF |
327 |
#endif /* SEAICE_ALLOW_DYNAMICS */ |
328 |
|
329 |
#ifdef ALLOW_AUTODIFF_TAMC |
330 |
CADJ STORE dwatn = comlev1, key=ikey_dynamics |
331 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
332 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
333 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
334 |
|
335 |
C Calculate ocean surface stress |
336 |
CALL SEAICE_OCEAN_STRESS ( myTime, myIter, myThid ) |
337 |
|
338 |
#ifdef SEAICE_ALLOW_DYNAMICS |
339 |
IF ( SEAICEuseDYNAMICS .AND. SEAICE_clipVelocities) THEN |
340 |
#ifdef ALLOW_AUTODIFF_TAMC |
341 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
342 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
343 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
344 |
c Put a cap on ice velocity |
345 |
c limit velocity to 0.40 m s-1 to avoid potential CFL violations |
346 |
c in open water areas (drift of zero thickness ice) |
347 |
DO bj=myByLo(myThid),myByHi(myThid) |
348 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
349 |
DO j=1-OLy,sNy+OLy |
350 |
DO i=1-OLx,sNx+OLx |
351 |
uIce(i,j,1,bi,bj)= |
352 |
& MAX(MIN(uIce(i,j,1,bi,bj),0.40 _d +00),-0.40 _d +00) |
353 |
vIce(i,j,1,bi,bj)= |
354 |
& MAX(MIN(vIce(i,j,1,bi,bj),0.40 _d +00),-0.40 _d +00) |
355 |
ENDDO |
356 |
ENDDO |
357 |
ENDDO |
358 |
ENDDO |
359 |
ENDIF |
360 |
#endif /* SEAICE_ALLOW_DYNAMICS */ |
361 |
|
362 |
#endif /* SEAICE_CGRID */ |
363 |
RETURN |
364 |
END |