1 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/dynsolver.F,v 1.44 2014/10/20 03:20:57 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "SEAICE_OPTIONS.h" |
5 |
#ifdef ALLOW_AUTODIFF |
6 |
# include "AUTODIFF_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
CStartOfInterface |
10 |
SUBROUTINE DYNSOLVER( myTime, myIter, myThid ) |
11 |
C *==========================================================* |
12 |
C | SUBROUTINE DYNSOLVER | |
13 |
C | o Ice dynamics using LSR solver | |
14 |
C | Zhang and Hibler, JGR, 102, 8691-8702, 1997 | |
15 |
C *==========================================================* |
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 "DYNVARS.h" |
25 |
#include "FFIELDS.h" |
26 |
#include "SEAICE_SIZE.h" |
27 |
#include "SEAICE_PARAMS.h" |
28 |
#include "SEAICE.h" |
29 |
#ifdef ALLOW_EXF |
30 |
# include "EXF_OPTIONS.h" |
31 |
# include "EXF_FIELDS.h" |
32 |
#endif |
33 |
#ifdef ALLOW_AUTODIFF_TAMC |
34 |
# include "tamc.h" |
35 |
#endif |
36 |
|
37 |
C === Routine arguments === |
38 |
C myTime - Simulation time |
39 |
C myIter - Simulation timestep number |
40 |
C myThid - Thread no. that called this routine. |
41 |
_RL myTime |
42 |
INTEGER myIter |
43 |
INTEGER myThid |
44 |
CEndOfInterface |
45 |
|
46 |
#ifndef SEAICE_CGRID |
47 |
|
48 |
#ifdef EXPLICIT_SSH_SLOPE |
49 |
#include "SURFACE.h" |
50 |
_RL phiSurf(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
51 |
#endif |
52 |
|
53 |
C === Local variables === |
54 |
C i,j,bi,bj - Loop counters |
55 |
|
56 |
INTEGER i, j, bi, bj |
57 |
_RL RHOICE, RHOAIR |
58 |
_RL COSWIN |
59 |
_RS SINWIN |
60 |
_RL ECCEN, ECM2, PSTAR, AAA |
61 |
_RL U1, V1 |
62 |
|
63 |
_RL COR_ICE (1-OLx:sNx+OLx,1-OLy:sNy+OLy, nSx,nSy) |
64 |
|
65 |
C-- FIRST SET UP BASIC CONSTANTS |
66 |
RHOICE=SEAICE_rhoIce |
67 |
RHOAIR=SEAICE_rhoAir |
68 |
ECCEN=SEAICE_eccen |
69 |
ECM2=ONE/(ECCEN**2) |
70 |
PSTAR=SEAICE_strength |
71 |
|
72 |
C-- introduce turning angle (default is zero) |
73 |
SINWIN=SIN(SEAICE_airTurnAngle*deg2rad) |
74 |
COSWIN=COS(SEAICE_airTurnAngle*deg2rad) |
75 |
|
76 |
C-- Compute proxy for geostrophic velocity, |
77 |
DO bj=myByLo(myThid),myByHi(myThid) |
78 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
79 |
DO j=0,sNy+1 |
80 |
DO i=0,sNx+1 |
81 |
GWATX(I,J,bi,bj)=HALF*(uVel(i,j,KGEO(I,J,bi,bj),bi,bj) |
82 |
& +uVel(i,j-1,KGEO(I,J,bi,bj),bi,bj)) |
83 |
GWATY(I,J,bi,bj)=HALF*(vVel(i,j,KGEO(I,J,bi,bj),bi,bj) |
84 |
& +vVel(i-1,j,KGEO(I,J,bi,bj),bi,bj)) |
85 |
#ifdef SEAICE_DEBUG |
86 |
c write(*,'(2i4,2i2,f7.1,7f12.3)') |
87 |
c & ,i,j,bi,bj,UVM(I,J,bi,bj) |
88 |
c & ,GWATX(I,J,bi,bj),GWATY(I,J,bi,bj) |
89 |
c & ,uVel(i+1,j,3,bi,bj),uVel(i+1,j+1,3,bi,bj) |
90 |
c & ,vVel(i,j+1,3,bi,bj),vVel(i+1,j+1,3,bi,bj) |
91 |
#endif |
92 |
ENDDO |
93 |
ENDDO |
94 |
ENDDO |
95 |
ENDDO |
96 |
|
97 |
C-- NOW SET UP MASS PER UNIT AREA AND CORIOLIS TERM |
98 |
DO bj=myByLo(myThid),myByHi(myThid) |
99 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
100 |
DO j=1-OLy,sNy+OLy |
101 |
DO i=1-OLx,sNx+OLx |
102 |
COR_ICE(i,j,bi,bj) = 0. |
103 |
ENDDO |
104 |
ENDDO |
105 |
DO j=1,sNy |
106 |
DO i=1,sNx |
107 |
AMASS(I,J,bi,bj)=RHOICE*QUART*( |
108 |
& HEFF(i,j ,bi,bj) + HEFF(i-1,j ,bi,bj) |
109 |
& +HEFF(i,j-1,bi,bj) + HEFF(i-1,j-1,bi,bj) ) |
110 |
COR_ICE(I,J,bi,bj)=AMASS(I,J,bi,bj) * _fCoriG(I,J,bi,bj) |
111 |
ENDDO |
112 |
ENDDO |
113 |
ENDDO |
114 |
ENDDO |
115 |
|
116 |
C-- NOW SET UP FORCING FIELDS |
117 |
|
118 |
C-- Wind stress is computed on South-West B-grid U/V |
119 |
C locations from wind on tracer locations |
120 |
DO bj=myByLo(myThid),myByHi(myThid) |
121 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
122 |
DO j=1,sNy |
123 |
DO i=1,sNx |
124 |
U1=QUART*(UWIND(I-1,J-1,bi,bj)+UWIND(I-1,J,bi,bj) |
125 |
& +UWIND(I ,J-1,bi,bj)+UWIND(I ,J,bi,bj)) |
126 |
V1=QUART*(VWIND(I-1,J-1,bi,bj)+VWIND(I-1,J,bi,bj) |
127 |
& +VWIND(I ,J-1,bi,bj)+VWIND(I ,J,bi,bj)) |
128 |
AAA=U1**2+V1**2 |
129 |
IF ( AAA .LE. SEAICE_EPS_SQ ) THEN |
130 |
AAA=SEAICE_EPS |
131 |
ELSE |
132 |
AAA=SQRT(AAA) |
133 |
ENDIF |
134 |
C first ocean surface stress |
135 |
DAIRN(I,J,bi,bj)=RHOAIR*OCEAN_drag |
136 |
& *(2.70 _d 0+0.142 _d 0*AAA+0.0764 _d 0*AAA*AAA) |
137 |
WINDX(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
138 |
& (COSWIN*U1-SIGN(SINWIN, _fCori(I,J,bi,bj))*V1) |
139 |
WINDY(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
140 |
& (SIGN(SINWIN, _fCori(I,J,bi,bj))*U1+COSWIN*V1) |
141 |
|
142 |
C now ice surface stress |
143 |
IF ( YC(I,J,bi,bj) .LT. ZERO ) THEN |
144 |
DAIRN(I,J,bi,bj) = |
145 |
& RHOAIR*(SEAICE_drag_south*AAA*AREA(I,J,bi,bj) |
146 |
& +OCEAN_drag*(2.70 _d 0+0.142 _d 0*AAA |
147 |
& +0.0764 _d 0*AAA*AAA)*(ONE-AREA(I,J,bi,bj))) |
148 |
ELSE |
149 |
DAIRN(I,J,bi,bj) = |
150 |
& RHOAIR*(SEAICE_drag*AAA*AREA(I,J,bi,bj) |
151 |
& +OCEAN_drag*(2.70 _d 0+0.142 _d 0*AAA |
152 |
& +0.0764 _d 0*AAA*AAA)*(ONE-AREA(I,J,bi,bj))) |
153 |
ENDIF |
154 |
FORCEX(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
155 |
& (COSWIN*U1-SIGN(SINWIN, _fCori(I,J,bi,bj))*V1) |
156 |
FORCEY(I,J,bi,bj)=DAIRN(I,J,bi,bj)* |
157 |
& (SIGN(SINWIN, _fCori(I,J,bi,bj))*U1+COSWIN*V1) |
158 |
ENDDO |
159 |
ENDDO |
160 |
ENDDO |
161 |
ENDDO |
162 |
|
163 |
DO bj=myByLo(myThid),myByHi(myThid) |
164 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
165 |
#ifdef EXPLICIT_SSH_SLOPE |
166 |
C-- Compute surface pressure at z==0: |
167 |
C- use actual sea surface height for tilt computations |
168 |
DO j=1-OLy,sNy+OLy |
169 |
DO i=1-OLx,sNx+OLx |
170 |
phiSurf(i,j) = Bo_surf(i,j,bi,bj)*etaN(i,j,bi,bj) |
171 |
ENDDO |
172 |
ENDDO |
173 |
#ifdef ATMOSPHERIC_LOADING |
174 |
C- add atmospheric loading and Sea-Ice loading |
175 |
IF ( useRealFreshWaterFlux ) THEN |
176 |
DO j=1-OLy,sNy+OLy |
177 |
DO i=1-OLx,sNx+OLx |
178 |
phiSurf(i,j) = phiSurf(i,j) |
179 |
& + ( pload(i,j,bi,bj) |
180 |
& +sIceLoad(i,j,bi,bj)*gravity |
181 |
& )*recip_rhoConst |
182 |
ENDDO |
183 |
ENDDO |
184 |
ELSE |
185 |
DO j=1-OLy,sNy+OLy |
186 |
DO i=1-OLx,sNx+OLx |
187 |
phiSurf(i,j) = phiSurf(i,j) |
188 |
& + pload(i,j,bi,bj)*recip_rhoConst |
189 |
ENDDO |
190 |
ENDDO |
191 |
ENDIF |
192 |
#endif /* ATMOSPHERIC_LOADING */ |
193 |
DO j=1-OLy+1,sNy+OLy |
194 |
DO i=1-OLx+1,sNx+OLx |
195 |
C-- NOW ADD IN TILT |
196 |
FORCEX(I,J,bi,bj)=FORCEX(I,J,bi,bj) |
197 |
& -AMASS(I,J,bi,bj) |
198 |
& *( (phiSurf(i, j )-phiSurf(i-1, j ))*maskW(i, j ,1,bi,bj) |
199 |
& +(phiSurf(i,j-1)-phiSurf(i-1,j-1))*maskW(i,j-1,1,bi,bj) |
200 |
& )*HALF*_recip_dxV(I,J,bi,bj) |
201 |
FORCEY(I,J,bi,bj)=FORCEY(I,J,bi,bj) |
202 |
& -AMASS(I,J,bi,bj) |
203 |
& *( (phiSurf( i ,j)-phiSurf( i ,j-1))*maskS( i ,j,1,bi,bj) |
204 |
& +(phiSurf(i-1,j)-phiSurf(i-1,j-1))*maskS(i-1,j,1,bi,bj) |
205 |
& )*HALF*_recip_dyU(I,J,bi,bj) |
206 |
C NOW KEEP FORCEX0 |
207 |
FORCEX0(I,J,bi,bj)=FORCEX(I,J,bi,bj) |
208 |
FORCEY0(I,J,bi,bj)=FORCEY(I,J,bi,bj) |
209 |
ENDDO |
210 |
ENDDO |
211 |
#endif /* EXPLICIT_SSH_SLOPE */ |
212 |
DO j=1-OLy,sNy+OLy |
213 |
DO i=1-OLx,sNx+OLx |
214 |
#ifndef EXPLICIT_SSH_SLOPE |
215 |
C-- NOW ADD IN TILT |
216 |
FORCEX(I,J,bi,bj)=FORCEX(I,J,bi,bj) |
217 |
& -COR_ICE(I,J,bi,bj)*GWATY(I,J,bi,bj) |
218 |
FORCEY(I,J,bi,bj)=FORCEY(I,J,bi,bj) |
219 |
& +COR_ICE(I,J,bi,bj)*GWATX(I,J,bi,bj) |
220 |
C NOW KEEP FORCEX0 |
221 |
FORCEX0(I,J,bi,bj)=FORCEX(I,J,bi,bj) |
222 |
FORCEY0(I,J,bi,bj)=FORCEY(I,J,bi,bj) |
223 |
#endif /* EXPLICIT_SSH_SLOPE */ |
224 |
C-- NOW SET UP ICE PRESSURE AND VISCOSITIES |
225 |
PRESS0(I,J,bi,bj)=PSTAR*HEFF(I,J,bi,bj) |
226 |
& *EXP(-SEAICE_cStar*(ONE-AREA(i,j,bi,bj))) |
227 |
CML ZMAX(I,J,bi,bj)=(5.0 _d +12/(2.0 _d +04))*PRESS0(I,J,bi,bj) |
228 |
ZMAX(I,J,bi,bj)=SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj) |
229 |
CML ZMIN(I,J,bi,bj)=4.0 _d +08 |
230 |
ZMIN(I,J,bi,bj)=SEAICE_zetaMin |
231 |
PRESS0(I,J,bi,bj)=PRESS0(I,J,bi,bj)*HEFFM(I,J,bi,bj) |
232 |
ENDDO |
233 |
ENDDO |
234 |
ENDDO |
235 |
ENDDO |
236 |
|
237 |
#ifdef SEAICE_ALLOW_DYNAMICS |
238 |
|
239 |
IF ( SEAICEuseDYNAMICS ) THEN |
240 |
|
241 |
#ifdef ALLOW_AUTODIFF_TAMC |
242 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
243 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
244 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
245 |
|
246 |
crg what about ETA,ZETA |
247 |
|
248 |
crg later c$taf loop = iteration uice,vice |
249 |
|
250 |
cdm c$taf store uice,vice = comlev1_seaice_ds, |
251 |
cdm c$taf& key = kii + (ikey_dynamics-1) |
252 |
C NOW DO PREDICTOR TIME STEP |
253 |
DO bj=myByLo(myThid),myByHi(myThid) |
254 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
255 |
DO j=1-OLy,sNy+OLy |
256 |
DO i=1-OLx,sNx+OLx |
257 |
UICENM1(I,J,bi,bj)=UICE(I,J,bi,bj) |
258 |
VICENM1(I,J,bi,bj)=VICE(I,J,bi,bj) |
259 |
UICEC(I,J,bi,bj)=UICE(I,J,bi,bj) |
260 |
VICEC(I,J,bi,bj)=VICE(I,J,bi,bj) |
261 |
ENDDO |
262 |
ENDDO |
263 |
ENDDO |
264 |
ENDDO |
265 |
|
266 |
C NOW LSR SCHEME (ZHANG-J/HIBLER 1997) |
267 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
268 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
269 |
CALL LSR( 1, myThid ) |
270 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
271 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
272 |
|
273 |
C NOW DO MODIFIED EULER STEP |
274 |
DO bj=myByLo(myThid),myByHi(myThid) |
275 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
276 |
DO j=1-OLy,sNy+OLy |
277 |
DO i=1-OLx,sNx+OLx |
278 |
UICE(I,J,bi,bj)=HALF*(UICE(I,J,bi,bj)+UICENM1(I,J,bi,bj)) |
279 |
VICE(I,J,bi,bj)=HALF*(VICE(I,J,bi,bj)+VICENM1(I,J,bi,bj)) |
280 |
UICEC(I,J,bi,bj)=UICE(I,J,bi,bj) |
281 |
VICEC(I,J,bi,bj)=VICE(I,J,bi,bj) |
282 |
ENDDO |
283 |
ENDDO |
284 |
ENDDO |
285 |
ENDDO |
286 |
|
287 |
C NOW LSR SCHEME (ZHANG-J/HIBLER 1997) |
288 |
CALL LSR( 2, myThid ) |
289 |
|
290 |
cdm c$taf store uice,vice = comlev1, key=ikey_dynamics |
291 |
|
292 |
ENDIF |
293 |
#endif /* SEAICE_ALLOW_DYNAMICS */ |
294 |
|
295 |
#ifdef ALLOW_DIAGNOSTICS |
296 |
IF ( useDiagnostics ) THEN |
297 |
CALL DIAGNOSTICS_FILL(zeta ,'SIzeta ',0,1 ,0,1,1,myThid) |
298 |
CALL DIAGNOSTICS_FILL(eta ,'SIeta ',0,1 ,0,1,1,myThid) |
299 |
CALL DIAGNOSTICS_FILL(press ,'SIpress ',0,1 ,0,1,1,myThid) |
300 |
ENDIF |
301 |
#endif /* ALLOW_DIAGNOSTICS */ |
302 |
|
303 |
C Calculate ocean surface stress |
304 |
CALL OSTRES ( COR_ICE, myThid ) |
305 |
|
306 |
#ifdef SEAICE_ALLOW_DYNAMICS |
307 |
#ifdef SEAICE_ALLOW_CLIPVELS |
308 |
IF ( SEAICEuseDYNAMICS .AND. SEAICE_clipVelocities) THEN |
309 |
#ifdef ALLOW_AUTODIFF_TAMC |
310 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
311 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
312 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
313 |
c Put a cap on ice velocity |
314 |
c limit velocity to 0.40 m s-1 to avoid potential CFL violations |
315 |
c in open water areas (drift of zero thickness ice) |
316 |
DO bj=myByLo(myThid),myByHi(myThid) |
317 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
318 |
DO j=1-OLy,sNy+OLy |
319 |
DO i=1-OLx,sNx+OLx |
320 |
#ifdef SEAICE_DEBUG |
321 |
c write(*,'(2i4,2i2,f7.1,7f12.3)') |
322 |
c & i,j,bi,bj,UVM(I,J,bi,bj),amass(i,j,bi,bj) |
323 |
c & ,gwatx(I,J,bi,bj),gwaty(i,j,bi,bj) |
324 |
c & ,forcex(I,J,bi,bj),forcey(i,j,bi,bj) |
325 |
c & ,uice(i,j,bi,bj) |
326 |
c & ,vice(i,j,bi,bj) |
327 |
#endif /* SEAICE_DEBUG */ |
328 |
UICE(i,j,bi,bj)=min(UICE(i,j,bi,bj),0.40 _d +00) |
329 |
VICE(i,j,bi,bj)=min(VICE(i,j,bi,bj),0.40 _d +00) |
330 |
ENDDO |
331 |
ENDDO |
332 |
ENDDO |
333 |
ENDDO |
334 |
#ifdef ALLOW_AUTODIFF_TAMC |
335 |
CADJ STORE uice = comlev1, key=ikey_dynamics |
336 |
CADJ STORE vice = comlev1, key=ikey_dynamics |
337 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
338 |
DO bj=myByLo(myThid),myByHi(myThid) |
339 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
340 |
DO j=1-OLy,sNy+OLy |
341 |
DO i=1-OLx,sNx+OLx |
342 |
UICE(i,j,bi,bj)=max(UICE(i,j,bi,bj),-0.40 _d +00) |
343 |
VICE(i,j,bi,bj)=max(VICE(i,j,bi,bj),-0.40 _d +00) |
344 |
ENDDO |
345 |
ENDDO |
346 |
ENDDO |
347 |
ENDDO |
348 |
|
349 |
ENDIF |
350 |
#endif /* SEAICE_ALLOW_CLIPVELS */ |
351 |
#endif /* SEAICE_ALLOW_DYNAMICS */ |
352 |
#endif /* NOT SEAICE_CGRID */ |
353 |
|
354 |
RETURN |
355 |
END |