/[MITgcm]/MITgcm/pkg/seaice/seaice_dynsolver.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_dynsolver.F

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


Revision 1.55 - (show annotations) (download)
Mon Apr 7 14:50:31 2014 UTC (10 years ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64w, checkpoint64v, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65d, checkpoint65e
Changes since 1.54: +3 -23 lines
add routines that calculate ice strength based on energy dissipation
during ridging following Rothrock (1975), Bitz et al (2001), Lipscomb
et al. (2007)

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_dynsolver.F,v 1.54 2012/10/21 04:34:07 heimbach Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: SEAICE_DYNSOLVER
8 C !INTERFACE:
9 SUBROUTINE SEAICE_DYNSOLVER( myTime, myIter, myThid )
10
11 C *==========================================================*
12 C | SUBROUTINE SEAICE_DYNSOLVER
13 C | o Ice dynamics using LSR solver
14 C | Zhang and Hibler, JGR, 102, 8691-8702, 1997
15 C | or EVP explicit solver by Hunke and Dukowicz, JPO 27,
16 C | 1849-1867 (1997)
17 C *==========================================================*
18 C | written by Martin Losch, March 2006
19 C *==========================================================*
20 C \ev
21
22 C !USES:
23 IMPLICIT NONE
24
25 C === Global variables ===
26 #include "SIZE.h"
27 #include "EEPARAMS.h"
28 #include "PARAMS.h"
29 #include "GRID.h"
30 #include "SURFACE.h"
31 #include "DYNVARS.h"
32 #include "FFIELDS.h"
33 #include "SEAICE_SIZE.h"
34 #include "SEAICE_PARAMS.h"
35 #include "SEAICE.h"
36
37 #ifdef ALLOW_AUTODIFF_TAMC
38 # include "tamc.h"
39 #endif
40
41 C !INPUT/OUTPUT PARAMETERS:
42 C === Routine arguments ===
43 C myTime :: Simulation time
44 C myIter :: Simulation timestep number
45 C myThid :: my Thread Id. number
46 _RL myTime
47 INTEGER myIter
48 INTEGER myThid
49 CEOP
50
51 #ifdef SEAICE_CGRID
52
53 C !FUNCTIONS:
54 LOGICAL DIFFERENT_MULTIPLE
55 EXTERNAL DIFFERENT_MULTIPLE
56
57 C !LOCAL VARIABLES:
58 C === Local variables ===
59 C i,j,bi,bj :: Loop counters
60 INTEGER i, j, bi, bj
61 _RL PSTAR
62 _RL phiSurf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63 _RL mask_uice, mask_vice
64
65 # ifdef ALLOW_AUTODIFF_TAMC
66 C Following re-initialisation breaks some "artificial" AD dependencies
67 C incured by IF (DIFFERENT_MULTIPLE ... statement
68 PSTAR = SEAICE_strength
69 DO bj=myByLo(myThid),myByHi(myThid)
70 DO bi=myBxLo(myThid),myBxHi(myThid)
71 DO j=1-OLy+1,sNy+OLy
72 DO i=1-OLx+1,sNx+OLx
73 PRESS0(i,j,bi,bj) = PSTAR*HEFF(i,j,bi,bj)
74 & *EXP(-20.0 _d 0*(ONE-AREA(i,j,bi,bj)))
75 ZMAX(I,J,bi,bj) = SEAICE_zetaMaxFac*PRESS0(I,J,bi,bj)
76 ZMIN(i,j,bi,bj) = SEAICE_zetaMin
77 PRESS0(i,j,bi,bj) = PRESS0(i,j,bi,bj)*HEFFM(i,j,bi,bj)
78 TAUX(i,j,bi,bj) = 0. _d 0
79 TAUY(i,j,bi,bj) = 0. _d 0
80 #ifdef SEAICE_ALLOW_FREEDRIFT
81 uice_fd(i,j,bi,bj)= 0. _d 0
82 vice_fd(i,j,bi,bj)= 0. _d 0
83 #endif
84 ENDDO
85 ENDDO
86 ENDDO
87 ENDDO
88 C
89 CADJ STORE uice = comlev1, key=ikey_dynamics, kind=isbyte
90 CADJ STORE vice = comlev1, key=ikey_dynamics, kind=isbyte
91 CADJ STORE uicenm1 = comlev1, key=ikey_dynamics, kind=isbyte
92 CADJ STORE vicenm1 = comlev1, key=ikey_dynamics, kind=isbyte
93 # endif /* ALLOW_AUTODIFF_TAMC */
94
95 IF (
96 & DIFFERENT_MULTIPLE(SEAICE_deltaTdyn,myTime,SEAICE_deltaTtherm)
97 & ) THEN
98
99 # ifdef ALLOW_AUTODIFF_TAMC
100 # ifdef SEAICE_ALLOW_EVP
101 CADJ STORE press0 = comlev1, key=ikey_dynamics, kind=isbyte
102 CADJ STORE zmax = comlev1, key=ikey_dynamics, kind=isbyte
103 # endif
104 # endif /* ALLOW_AUTODIFF_TAMC */
105
106 C-- FIRST SET UP BASIC CONSTANTS
107 PSTAR = SEAICE_strength
108
109 C-- NOW SET UP MASS PER UNIT AREA AND CORIOLIS TERM
110 DO bj=myByLo(myThid),myByHi(myThid)
111 DO bi=myBxLo(myThid),myBxHi(myThid)
112 DO j=1-OLy+1,sNy+OLy
113 DO i=1-OLx+1,sNx+OLx
114 seaiceMassC(I,J,bi,bj)=SEAICE_rhoIce*HEFF(i,j,bi,bj)
115 seaiceMassU(I,J,bi,bj)=SEAICE_rhoIce*HALF*(
116 & HEFF(i,j,bi,bj) + HEFF(i-1,j ,bi,bj) )
117 seaiceMassV(I,J,bi,bj)=SEAICE_rhoIce*HALF*(
118 & HEFF(i,j,bi,bj) + HEFF(i ,j-1,bi,bj) )
119 ENDDO
120 ENDDO
121 ENDDO
122 ENDDO
123
124 #ifndef ALLOW_AUTODIFF_TAMC
125 IF ( SEAICE_maskRHS ) THEN
126 C dynamic masking of areas with no ice
127 DO bj=myByLo(myThid),myByHi(myThid)
128 DO bi=myBxLo(myThid),myBxHi(myThid)
129 DO j=1-OLy+1,sNy+OLy
130 DO i=1-OLx+1,sNx+OLx
131 seaiceMaskU(I,J,bi,bj)=AREA(i,j,bi,bj)+AREA(I-1,J,bi,bj)
132 mask_uice=HEFFM(i,j,bi,bj)+HEFFM(i-1,j ,bi,bj)
133 IF ( (seaiceMaskU(I,J,bi,bj) .GT. 0. _d 0) .AND.
134 & (mask_uice .GT. 1.5 _d 0) ) THEN
135 seaiceMaskU(I,J,bi,bj) = 1. _d 0
136 ELSE
137 seaiceMaskU(I,J,bi,bj) = 0. _d 0
138 ENDIF
139 seaiceMaskV(I,J,bi,bj)=AREA(i,j,bi,bj)+AREA(I,J-1,bi,bj)
140 mask_vice=HEFFM(i,j,bi,bj)+HEFFM(i ,j-1,bi,bj)
141 IF ( (seaiceMaskV(I,J,bi,bj) .GT. 0. _d 0) .AND.
142 & (mask_vice .GT. 1.5 _d 0) ) THEN
143 seaiceMaskV(I,J,bi,bj) = 1. _d 0
144 ELSE
145 seaiceMaskV(I,J,bi,bj) = 0. _d 0
146 ENDIF
147 ENDDO
148 ENDDO
149 ENDDO
150 ENDDO
151 CALL EXCH_UV_XY_RL( seaiceMaskU, seaiceMaskV, .FALSE., myThid )
152 ENDIF
153 #endif /* ndef ALLOW_AUTODIFF_TAMC */
154
155 C-- NOW SET UP FORCING FIELDS
156
157 C initialise fields
158 DO bj=myByLo(myThid),myByHi(myThid)
159 DO bi=myBxLo(myThid),myBxHi(myThid)
160 DO j=1-OLy,sNy+OLy
161 DO i=1-OLx,sNx+OLx
162 TAUX (I,J,bi,bj)= 0. _d 0
163 TAUY (I,J,bi,bj)= 0. _d 0
164 #ifdef ALLOW_AUTODIFF_TAMC
165 # ifdef SEAICE_ALLOW_EVP
166 stressDivergenceX(I,J,bi,bj) = 0. _d 0
167 stressDivergenceY(I,J,bi,bj) = 0. _d 0
168 # endif
169 #endif
170 ENDDO
171 ENDDO
172 ENDDO
173 ENDDO
174
175 # ifdef ALLOW_AUTODIFF_TAMC
176 CADJ STORE uice = comlev1, key=ikey_dynamics, kind=isbyte
177 CADJ STORE vice = comlev1, key=ikey_dynamics, kind=isbyte
178 # endif /* ALLOW_AUTODIFF_TAMC */
179 C-- interface of dynamics with atmopheric forcing fields (wind/stress)
180 CALL SEAICE_GET_DYNFORCING (
181 I uIce, vIce,
182 O TAUX, TAUY,
183 I myTime, myIter, myThid )
184
185 DO bj=myByLo(myThid),myByHi(myThid)
186 DO bi=myBxLo(myThid),myBxHi(myThid)
187 C-- Compute surface pressure at z==0:
188 C- use actual sea surface height for tilt computations
189 DO j=1-OLy,sNy+OLy
190 DO i=1-OLx,sNx+OLx
191 phiSurf(i,j) = Bo_surf(i,j,bi,bj)*etaN(i,j,bi,bj)
192 ENDDO
193 ENDDO
194 #ifdef ATMOSPHERIC_LOADING
195 C- add atmospheric loading and Sea-Ice loading
196 IF ( useRealFreshWaterFlux ) THEN
197 DO j=1-OLy,sNy+OLy
198 DO i=1-OLx,sNx+OLx
199 phiSurf(i,j) = phiSurf(i,j)
200 & + ( pload(i,j,bi,bj)
201 & +sIceLoad(i,j,bi,bj)*gravity
202 & )*recip_rhoConst
203 ENDDO
204 ENDDO
205 ELSE
206 DO j=1-OLy,sNy+OLy
207 DO i=1-OLx,sNx+OLx
208 phiSurf(i,j) = phiSurf(i,j)
209 & + pload(i,j,bi,bj)*recip_rhoConst
210 ENDDO
211 ENDDO
212 ENDIF
213 #endif /* ATMOSPHERIC_LOADING */
214 DO j=1-OLy+1,sNy+OLy
215 DO i=1-OLx+1,sNx+OLx
216 C-- basic forcing by wind stress
217 FORCEX0(I,J,bi,bj)=TAUX(I,J,bi,bj)
218 FORCEY0(I,J,bi,bj)=TAUY(I,J,bi,bj)
219 ENDDO
220 ENDDO
221
222 IF ( SEAICEuseTILT ) then
223 DO j=1-OLy+1,sNy+OLy
224 DO i=1-OLx+1,sNx+OLx
225 C-- now add in tilt
226 FORCEX0(I,J,bi,bj)=FORCEX0(I,J,bi,bj)
227 & -seaiceMassU(I,J,bi,bj)*_recip_dxC(I,J,bi,bj)
228 & *( phiSurf(i,j)-phiSurf(i-1,j) )
229 FORCEY0(I,J,bi,bj)=FORCEY0(I,J,bi,bj)
230 & -seaiceMassV(I,J,bi,bj)* _recip_dyC(I,J,bi,bj)
231 & *( phiSurf(i,j)-phiSurf(i,j-1) )
232 ENDDO
233 ENDDO
234 ENDIF
235
236 CALL SEAICE_CALC_ICE_STRENGTH( bi, bj, myTime, myIter, myThid )
237
238 ENDDO
239 ENDDO
240
241 #ifdef SEAICE_ALLOW_DYNAMICS
242 IF ( SEAICEuseDYNAMICS ) THEN
243
244 #ifdef SEAICE_ALLOW_FREEDRIFT
245 IF ( SEAICEuseFREEDRIFT .OR. SEAICEuseEVP
246 & .OR. LSR_mixIniGuess.EQ.0 ) THEN
247 CALL SEAICE_FREEDRIFT( myTime, myIter, myThid )
248 ENDIF
249 IF ( SEAICEuseFREEDRIFT ) THEN
250 DO bj=myByLo(myThid),myByHi(myThid)
251 DO bi=myBxLo(myThid),myBxHi(myThid)
252 DO j=1-OLy,sNy+OLy
253 DO i=1-OLx,sNx+OLx
254 uIce(i,j,bi,bj) = uIce_fd(i,j,bi,bj)
255 vIce(i,j,bi,bj) = vIce_fd(i,j,bi,bj)
256 stressDivergenceX(i,j,bi,bj) = 0. _d 0
257 stressDivergenceY(i,j,bi,bj) = 0. _d 0
258 ENDDO
259 ENDDO
260 ENDDO
261 ENDDO
262 ENDIF
263 #endif /* SEAICE_ALLOW_FREEDRIFT */
264
265 #ifdef ALLOW_OBCS
266 IF ( useOBCS ) THEN
267 CALL OBCS_APPLY_UVICE( uIce, vIce, myThid )
268 ENDIF
269 #endif /* ALLOW_OBCS */
270
271 #ifdef SEAICE_ALLOW_EVP
272 # ifdef ALLOW_AUTODIFF_TAMC
273 CADJ STORE uice = comlev1, key=ikey_dynamics, kind=isbyte
274 CADJ STORE vice = comlev1, key=ikey_dynamics, kind=isbyte
275 CADJ STORE uicenm1 = comlev1, key=ikey_dynamics, kind=isbyte
276 CADJ STORE vicenm1 = comlev1, key=ikey_dynamics, kind=isbyte
277 # endif /* ALLOW_AUTODIFF_TAMC */
278 IF ( SEAICEuseEVP ) THEN
279 C Elastic-Viscous-Plastic solver, following Hunke (2001)
280 CALL SEAICE_EVP( myTime, myIter, myThid )
281 ENDIF
282 #endif /* SEAICE_ALLOW_EVP */
283
284 #ifdef SEAICE_ALLOW_JFNK
285 IF ( SEAICEuseJFNK ) THEN
286 # ifdef ALLOW_AUTODIFF_TAMC
287 STOP 'Adjoint does not work with JFNK solver.'
288 # else
289 C Jacobian-free Newton Krylov solver (Lemieux et al. 2010, 2012)
290 CALL SEAICE_JFNK( myTime, myIter, myThid )
291 # endif /* ALLOW_AUTODIFF_TAMC */
292 ENDIF
293 #endif /* SEAICE_ALLOW_JFNK */
294
295 #if defined(SEAICE_ALLOW_EVP) || defined(SEAICE_ALLOW_JFNK) \
296 || defined(SEAICE_ALLOW_FREEDRIFT)
297 IF ( .NOT.SEAICEuseFREEDRIFT .AND. .NOT.SEAICEuseEVP
298 & .AND. .NOT.SEAICEuseJFNK ) THEN
299 #endif /* SEAICE_ALLOW_EVP or SEAICE_ALLOW_FREEDRIFT */
300 C LSR scheme (Zhang-J/Hibler 1997), ported to a C-grid
301 CALL SEAICE_LSR( myTime, myIter, myThid )
302 #if defined(SEAICE_ALLOW_EVP) || defined(SEAICE_ALLOW_JFNK) \
303 || defined(SEAICE_ALLOW_FREEDRIFT)
304 ENDIF
305 #endif /* SEAICE_ALLOW_EVP or SEAICE_ALLOW_FREEDRIFT */
306
307 C End of IF (SEAICEuseDYNAMICS ...
308 ENDIF
309 #endif /* SEAICE_ALLOW_DYNAMICS */
310
311 C End of IF (DIFFERENT_MULTIPLE ...
312 ENDIF
313
314 #ifdef ALLOW_AUTODIFF_TAMC
315 CADJ STORE uice = comlev1, key=ikey_dynamics, kind=isbyte
316 CADJ STORE vice = comlev1, key=ikey_dynamics, kind=isbyte
317 CADJ STORE stressDivergenceX = comlev1,
318 CADJ & key=ikey_dynamics, kind=isbyte
319 CADJ STORE stressDivergenceY = comlev1,
320 CADJ & key=ikey_dynamics, kind=isbyte
321 CADJ STORE DWATN = comlev1, key=ikey_dynamics, kind=isbyte
322 #endif /* ALLOW_AUTODIFF_TAMC */
323
324 C Calculate ocean surface stress
325 CALL SEAICE_OCEAN_STRESS ( myTime, myIter, myThid )
326
327 #ifdef SEAICE_ALLOW_DYNAMICS
328 #ifdef SEAICE_ALLOW_CLIPVELS
329 IF ( SEAICEuseDYNAMICS .AND. SEAICE_clipVelocities) THEN
330 #ifdef ALLOW_AUTODIFF_TAMC
331 CADJ STORE uice = comlev1, key=ikey_dynamics, kind=isbyte
332 CADJ STORE vice = comlev1, key=ikey_dynamics, kind=isbyte
333 #endif /* ALLOW_AUTODIFF_TAMC */
334 c Put a cap on ice velocity
335 c limit velocity to 0.40 m s-1 to avoid potential CFL violations
336 c in open water areas (drift of zero thickness ice)
337 DO bj=myByLo(myThid),myByHi(myThid)
338 DO bi=myBxLo(myThid),myBxHi(myThid)
339 DO j=1-OLy,sNy+OLy
340 DO i=1-OLx,sNx+OLx
341 uIce(i,j,bi,bj)=
342 & MAX(MIN(uIce(i,j,bi,bj),0.40 _d +00),-0.40 _d +00)
343 vIce(i,j,bi,bj)=
344 & MAX(MIN(vIce(i,j,bi,bj),0.40 _d +00),-0.40 _d +00)
345 ENDDO
346 ENDDO
347 ENDDO
348 ENDDO
349 ENDIF
350 #endif /* SEAICE_ALLOW_CLIPVELS */
351 #endif /* SEAICE_ALLOW_DYNAMICS */
352
353 #endif /* SEAICE_CGRID */
354 RETURN
355 END

  ViewVC Help
Powered by ViewVC 1.1.22