/[MITgcm]/MITgcm/pkg/thsice/thsice_calc_thickn.F
ViewVC logotype

Contents of /MITgcm/pkg/thsice/thsice_calc_thickn.F

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


Revision 1.15 - (show annotations) (download)
Fri Apr 11 21:31:14 2008 UTC (16 years, 2 months ago) by dfer
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59r, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q
Changes since 1.14: +16 -11 lines
- Now hsMax limits snow heitght by turning it into ice following
  flooding scheme (and now conserving energy)
- Slight rewriting of the formula in thsice_calc_thick.F
- Parameter rhowi eliminated to be replaced by floodFac = (rhosw-rhoi)/rhos
--> changes results of global_ocean.cs32x15.icedyn and global_ocean.cs32x15.thsice

1 C $Header: /u/gcmpack/MITgcm/pkg/thsice/thsice_calc_thickn.F,v 1.14 2007/09/17 21:58:15 jmc Exp $
2 C $Name: $
3
4 #include "THSICE_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: THSICE_CALC_THICKN
8 C !INTERFACE:
9 SUBROUTINE THSICE_CALC_THICKN(
10 I bi, bj, siLo, siHi, sjLo, sjHi,
11 I iMin,iMax, jMin,jMax, dBugFlag,
12 I iceMask, tFrz, tOce, v2oc,
13 I snowP, prcAtm, sHeat, flxCnB,
14 U icFrac, hIce, hSnow, tSrf, qIc1, qIc2,
15 U frwAtm, fzMlOc, flx2oc,
16 O frw2oc, fsalt,
17 I myTime, myIter, myThid )
18 C !DESCRIPTION: \bv
19 C *==========================================================*
20 C | S/R THSICE_CALC_THICKN
21 C | o Calculate ice & snow thickness changes
22 C *==========================================================*
23 C \ev
24 C ADAPTED FROM:
25 C LANL CICE.v2.0.2
26 C-----------------------------------------------------------------------
27 C.. thermodynamics (vertical physics) based on M. Winton 3-layer model
28 C.. See Bitz, C. M. and W. H. Lipscomb, 1999: "An energy-conserving
29 C.. thermodynamic sea ice model for climate study." J. Geophys.
30 C.. Res., 104, 15669 - 15677.
31 C.. Winton, M., 1999: "A reformulated three-layer sea ice model."
32 C.. Submitted to J. Atmos. Ocean. Technol.
33 C.. authors Elizabeth C. Hunke and William Lipscomb
34 C.. Fluid Dynamics Group, Los Alamos National Laboratory
35 C-----------------------------------------------------------------------
36 Cc****subroutine thermo_winton(n,fice,fsnow,dqice,dTsfc)
37 C.. Compute temperature change using Winton model with 2 ice layers, of
38 C.. which only the top layer has a variable heat capacity.
39
40 C---------------------------------
41 C parameters that control the partitioning between lateral (ice area) and
42 C vertical (ice thickness) ice volume changes.
43 C a) surface melting and bottom melting (runtime parameter: fracEnMelt):
44 C frace is the fraction of available heat that is used for
45 C lateral melting (and 1-frace reduces the thickness ) when
46 C o hi < hThinIce & frac > lowIcFrac2 : frace=1 (lateral melting only)
47 C o hThinIce < hi < hThickIce & frac > lowIcFrac1 : frace=fracEnMelt
48 C o hi > hThickIce or frac < lowIcFrac1 : frace=0 (thinning only)
49 C b) ocean freezing (and ice forming):
50 C - conductive heat flux (below sea-ice) always increases thickness.
51 C - under sea-ice, freezing potential (x iceFraction) is used to increase ice
52 C thickness or ice fraction (lateral growth), according to:
53 C o hi < hThinIce : use freezing potential to grow ice vertically;
54 C o hThinIce < hi < hThickIce : use partition factor fracEnFreez for lateral growth
55 c and (1-fracEnFreez) to increase thickness.
56 C o hi > hThickIce : use all freezing potential to grow ice laterally
57 C (up to areaMax)
58 C - over open ocean, use freezing potential [x(1-iceFraction)] to grow ice laterally
59 C - lateral growth forms ice of the same or =hNewIceMax thickness, the less of the 2.
60 C - starts to form sea-ice over fraction iceMaskMin, as minimum ice-volume is reached
61 C---------------------------------
62 C !USES:
63 IMPLICIT NONE
64
65 C == Global variables ===
66 #include "EEPARAMS.h"
67 #include "THSICE_SIZE.h"
68 #include "THSICE_PARAMS.h"
69 #ifdef ALLOW_AUTODIFF_TAMC
70 # include "SIZE.h"
71 # include "tamc.h"
72 # include "tamc_keys.h"
73 #endif
74
75 C !INPUT/OUTPUT PARAMETERS:
76 C == Routine Arguments ==
77 C siLo,siHi :: size of input/output array: 1rst dim. lower,higher bounds
78 C sjLo,sjHi :: size of input/output array: 2nd dim. lower,higher bounds
79 C bi,bj :: tile indices
80 C iMin,iMax :: computation domain: 1rst index range
81 C jMin,jMax :: computation domain: 2nd index range
82 C dBugFlag :: allow to print debugging stuff (e.g. on 1 grid point).
83 C--- Input:
84 C iceMask :: sea-ice fractional mask [0-1]
85 C tFrz (Tf) :: sea-water freezing temperature [oC] (function of S)
86 C tOce (oceTs) :: surface level oceanic temperature [oC]
87 C v2oc (oceV2s) :: square of ocean surface-level velocity [m2/s2]
88 C snowP (snowPr) :: snow precipitation [kg/m2/s]
89 C prcAtm (=) :: total precip from the atmosphere [kg/m2/s]
90 C sHeat(sHeating):: surf heating flux left to melt snow or ice (= Atmos-conduction)
91 C flxCnB (=) :: heat flux conducted through the ice to bottom surface
92 C--- Modified (input&output):
93 C icFrac(iceFrac):: fraction of grid area covered in ice
94 C hIce (hi) :: ice height [m]
95 C hSnow (hs) :: snow height [m]
96 C tSrf (Tsf) :: surface (ice or snow) temperature
97 C qIc1 (qicen) :: ice enthalpy (J/kg), 1rst level
98 C qIc2 (qicen) :: ice enthalpy (J/kg), 2nd level
99 C frwAtm (evpAtm):: evaporation to the atmosphere [kg/m2/s] (>0 if evaporate)
100 C fzMlOc (frzmlt):: ocean mixed-layer freezing/melting potential [W/m2]
101 C flx2oc (qleft) :: net heat flux to ocean [W/m2] (> 0 downward)
102 C--- Output
103 C frw2oc (fresh) :: Total fresh water flux to ocean [kg/m2/s] (> 0 downward)
104 C fsalt (=) :: salt flux to ocean [g/m2/s] (> 0 downward)
105 C--- Input:
106 C myTime :: current Time of simulation [s]
107 C myIter :: current Iteration number in simulation
108 C myThid :: my Thread Id number
109 INTEGER siLo, siHi, sjLo, sjHi
110 INTEGER bi,bj
111 INTEGER iMin, iMax
112 INTEGER jMin, jMax
113 LOGICAL dBugFlag
114 _RL iceMask(siLo:siHi,sjLo:sjHi)
115 _RL tFrz (siLo:siHi,sjLo:sjHi)
116 _RL tOce (siLo:siHi,sjLo:sjHi)
117 _RL v2oc (siLo:siHi,sjLo:sjHi)
118 _RL snowP (siLo:siHi,sjLo:sjHi)
119 _RL prcAtm (siLo:siHi,sjLo:sjHi)
120 _RL sHeat (siLo:siHi,sjLo:sjHi)
121 _RL flxCnB (siLo:siHi,sjLo:sjHi)
122 _RL icFrac (siLo:siHi,sjLo:sjHi)
123 _RL hIce (siLo:siHi,sjLo:sjHi)
124 _RL hSnow (siLo:siHi,sjLo:sjHi)
125 _RL tSrf (siLo:siHi,sjLo:sjHi)
126 _RL qIc1 (siLo:siHi,sjLo:sjHi)
127 _RL qIc2 (siLo:siHi,sjLo:sjHi)
128 _RL frwAtm (siLo:siHi,sjLo:sjHi)
129 _RL fzMlOc (siLo:siHi,sjLo:sjHi)
130 _RL flx2oc (siLo:siHi,sjLo:sjHi)
131 _RL frw2oc (siLo:siHi,sjLo:sjHi)
132 _RL fsalt (siLo:siHi,sjLo:sjHi)
133 _RL myTime
134 INTEGER myIter
135 INTEGER myThid
136 CEOP
137
138 #ifdef ALLOW_THSICE
139
140 C !LOCAL VARIABLES:
141 C--- local copy of input/output argument list variables (see description above)
142 _RL frzmlt
143 _RL Tf
144 _RL oceTs, oceV2s, snowPr
145 _RL sHeating
146 c _RL flxCnB
147 c _RL evpAtm
148 _RL iceFrac
149 _RL hi
150 _RL hs
151 _RL Tsf
152 _RL qicen(nlyr)
153 _RL qleft
154 _RL fresh
155 c _RL fsalt
156 C == Local Variables ==
157 INTEGER i,j,k ! loop indices
158 _RL rec_nlyr ! reciprocal of number of ice layers (real value)
159 C evap :: evaporation over snow/ice [kg/m2/s] (>0 if evaporate)
160 C Fbot :: oceanic heat flux used to melt/form ice [W/m2]
161 _RL evap
162 _RL Fbot
163 _RL etop ! energy for top melting (J m-2)
164 _RL ebot ! energy for bottom melting (J m-2)
165 _RL etope ! energy (from top) for lateral melting (J m-2)
166 _RL ebote ! energy (from bottom) for lateral melting (J m-2)
167 _RL extend ! total energy for lateral melting (J m-2)
168 _RL hnew(nlyr) ! new ice layer thickness (m)
169 _RL hlyr ! individual ice layer thickness (m)
170 _RL dhi ! change in ice thickness
171 _RL dhs ! change in snow thickness
172 _RL rq ! rho * q for a layer
173 _RL rqh ! rho * q * h for a layer
174 _RL qbot ! enthalpy for new ice at bottom surf (J/kg)
175 _RL dt ! timestep
176 _RL esurp ! surplus energy from melting (J m-2)
177 _RL mwater0 ! fresh water mass gained/lost (kg/m^2)
178 _RL msalt0 ! salt gained/lost (kg/m^2)
179 _RL freshe ! fresh water gain from extension melting
180 _RL salte ! salt gained from extension melting
181
182 _RL ustar, cpchr
183 _RL chi, chs
184 _RL frace, rs, hq
185 C lowIcFrac1 :: ice-fraction lower limit above which partial (lowIcFrac1)
186 C lowIcFrac2 :: or full (lowIcFrac2) lateral melting is allowed.
187 _RL lowIcFrac1, lowIcFrac2
188 #ifdef THSICE_FRACEN_POWERLAW
189 INTEGER powerLaw
190 _RL rec_pLaw
191 _RL c1Mlt, c2Mlt, aMlt, hMlt
192 _RL c1Frz, c2Frz, aFrz, hFrz
193 _RL enFrcMlt, xxMlt, tmpMlt
194 _RL enFrcFrz, xxFrz, tmpFrz
195 #endif
196
197 C- define grid-point location where to print debugging values
198 #include "THSICE_DEBUG.h"
199
200 1010 FORMAT(A,I3,3F8.3)
201 1020 FORMAT(A,1P4E11.3)
202
203 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
204
205 #ifdef ALLOW_AUTODIFF_TAMC
206 act1 = bi - myBxLo(myThid)
207 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
208 act2 = bj - myByLo(myThid)
209 max2 = myByHi(myThid) - myByLo(myThid) + 1
210 act3 = myThid - 1
211 max3 = nTx*nTy
212 act4 = ikey_dynamics - 1
213 #endif /* ALLOW_AUTODIFF_TAMC */
214
215 rec_nlyr = nlyr
216 rec_nlyr = 1. _d 0 / rec_nlyr
217 dt = thSIce_deltaT
218
219 C for now, use hard coded threshold (iceMaskMin +1.% and +10.%)
220 lowIcFrac1 = iceMaskMin*1.01 _d 0
221 lowIcFrac2 = iceMaskMin*1.10 _d 0
222 #ifdef THSICE_FRACEN_POWERLAW
223 IF ( powerLawExp2 .GE. 1 ) THEN
224 powerLaw = 1 + 2**powerLawExp2
225 rec_pLaw = powerLaw
226 rec_pLaw = 1. _d 0 / rec_pLaw
227 C- Coef for melting:
228 C lateral-melting energy fraction = fracEnMelt - [ aMlt*(hi-hMlt) ]^powerLaw
229 c1Mlt = fracEnMelt**rec_pLaw
230 c2Mlt = (1. _d 0 - fracEnMelt)**rec_pLaw
231 aMlt = (c1Mlt+c2Mlt)/(hThickIce-hThinIce)
232 hMlt = hThinIce+c2Mlt/aMlt
233 C- Coef for freezing:
234 C thickening energy fraction = fracEnFreez - [ aFrz*(hi-hFrz) ]^powerLaw
235 c1Frz = fracEnFreez**rec_pLaw
236 c2Frz = (1. _d 0 -fracEnFreez)**rec_pLaw
237 aFrz = (c1Frz+c2Frz)/(hThickIce-hThinIce)
238 hFrz = hThinIce+c2Frz/aFrz
239 ELSE
240 C- Linear relation
241 powerLaw = 1
242 aMlt = -1. _d 0 /(hThickIce-hThinIce)
243 hMlt = hThickIce
244 aFrz = -1. _d 0 /(hThickIce-hThinIce)
245 hFrz = hThickIce
246 ENDIF
247 #endif /* THSICE_FRACEN_POWERLAW */
248
249
250 DO j = jMin, jMax
251 DO i = iMin, iMax
252 #ifdef ALLOW_AUTODIFF_TAMC
253 ikey_1 = i
254 & + sNx*(j-1)
255 & + sNx*sNy*act1
256 & + sNx*sNy*max1*act2
257 & + sNx*sNy*max1*max2*act3
258 & + sNx*sNy*max1*max2*max3*act4
259 #endif /* ALLOW_AUTODIFF_TAMC */
260 C--
261 #ifdef ALLOW_AUTODIFF_TAMC
262 CADJ STORE frwatm(i,j) = comlev1_thsice_1, key=ikey_1
263 CADJ STORE fzmloc(i,j) = comlev1_thsice_1, key=ikey_1
264 CADJ STORE hice(i,j) = comlev1_thsice_1, key=ikey_1
265 CADJ STORE hsnow(i,j) = comlev1_thsice_1, key=ikey_1
266 CADJ STORE icfrac(i,j) = comlev1_thsice_1, key=ikey_1
267 CADJ STORE qic1(i,j) = comlev1_thsice_1, key=ikey_1
268 CADJ STORE qic2(i,j) = comlev1_thsice_1, key=ikey_1
269 #endif
270
271 IF (iceMask(i,j).GT.0. _d 0) THEN
272 Tf = tFrz(i,j)
273 oceTs = tOce(i,j)
274 oceV2s = v2oc(i,j)
275 snowPr = snowP(i,j)
276 c prcAtm = prcAtm(i,j)
277 sHeating= sHeat(i,j)
278 c flxCnB = flxCnB(i,j)
279 iceFrac = icFrac(i,j)
280 hi = hIce(i,j)
281 hs = hSnow(i,j)
282 Tsf = tSrf(i,j)
283 qicen(1)= qIc1(i,j)
284 qicen(2)= qIc2(i,j)
285 c evpAtm = frwAtm(i,j)
286 frzmlt = fzMlOc(i,j)
287 qleft = flx2oc(i,j)
288 c fresh = frw2oc(i,j)
289 c fsalt = fsalt(i,j)
290 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
291 C initialize energies
292 esurp = 0. _d 0
293
294 evap = frwAtm(i,j)
295
296 C......................................................................
297 C.. Compute growth and/or melting at the top and bottom surfaces.......
298 C......................................................................
299
300 #ifdef THSICE_FRACEN_POWERLAW
301 xxMlt = aMlt*(hi-hMlt)
302 xxFrz = aFrz*(hi-hFrz)
303 c--
304 IF ( powerLawExp2 .GE. 1 ) THEN
305 tmpMlt = xxMlt
306 tmpFrz = xxFrz
307 DO k=1,powerLawExp2
308 tmpMlt = tmpMlt*tmpMlt
309 tmpFrz = tmpFrz*tmpFrz
310 ENDDO
311 xxMlt = xxMlt*tmpMlt
312 xxFrz = xxFrz*tmpFrz
313 c xxMlt = xxMlt**powerLaw
314 c xxFrz = xxFrz**powerLaw
315 xxMlt = fracEnMelt -xxMlt
316 xxFrz = fracEnFreez-xxFrz
317 ENDIF
318 enFrcMlt = MAX( 0. _d 0, MIN( xxMlt, 1. _d 0 ) )
319 enFrcFrz = MAX( 0. _d 0, MIN( xxFrz, 1. _d 0 ) )
320 #endif /* THSICE_FRACEN_POWERLAW */
321
322 IF (frzmlt.GE. 0. _d 0) THEN
323 C !-----------------------------------------------------------------
324 C ! freezing conditions
325 C !-----------------------------------------------------------------
326 Fbot = frzmlt
327 IF ( iceFrac.LT.iceMaskMax ) THEN
328 #ifdef THSICE_FRACEN_POWERLAW
329 Fbot = enFrcFrz*frzmlt
330 #else /* THSICE_FRACEN_POWERLAW */
331 IF (hi.GT.hThickIce) THEN
332 C if higher than hThickIce, use all frzmlt energy to grow extra ice
333 Fbot = 0. _d 0
334 ELSEIF (hi.GE.hThinIce) THEN
335 C between hThinIce & hThickIce, use partition factor fracEnFreez
336 Fbot = (1. _d 0 - fracEnFreez)*frzmlt
337 ENDIF
338 #endif /* THSICE_FRACEN_POWERLAW */
339 ENDIF
340 ELSE
341 C !-----------------------------------------------------------------
342 C ! melting conditions
343 C !-----------------------------------------------------------------
344 ustar = 5. _d -2 !for no currents
345 C frictional velocity between ice and water
346 ustar = SQRT(0.00536 _d 0*oceV2s)
347 ustar=max(5. _d -3,ustar)
348 cpchr =cpWater*rhosw*bMeltCoef
349 Fbot = cpchr*(Tf-oceTs)*ustar ! < 0
350 Fbot = max(Fbot,frzmlt) ! frzmlt < Fbot < 0
351 Fbot = min(Fbot,0. _d 0)
352 ENDIF
353
354 C mass of fresh water and salt initially present in ice
355 mwater0 = rhos*hs + rhoi*hi
356 msalt0 = rhoi*hi*saltIce
357
358 #ifdef ALLOW_DBUG_THSICE
359 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
360 & 'ThSI_CALC_TH: evpAtm, frzmlt, Fbot =', frwAtm(i,j),frzmlt,Fbot
361 #endif
362
363 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
364
365 C Compute energy available for melting/growth.
366
367 #ifdef THSICE_FRACEN_POWERLAW
368 IF ( fracEnMelt.EQ.0. _d 0 ) THEN
369 frace = 0. _d 0
370 ELSE
371 frace = (iceFrac - lowIcFrac1)/(lowIcFrac2-iceMaskMin)
372 frace = MIN( enFrcMlt, MAX( 0. _d 0, frace ) )
373 ENDIF
374 #else /* THSICE_FRACEN_POWERLAW */
375 IF ( hi.GT.hThickIce .OR. fracEnMelt.EQ.0. _d 0 ) THEN
376 C above certain height (or when no ice fractionation), only melt from top
377 frace = 0. _d 0
378 ELSEIF (hi.LT.hThinIce) THEN
379 C below a certain height, all energy goes to changing ice extent
380 frace = 1. _d 0
381 ELSE
382 frace = fracEnMelt
383 ENDIF
384 C Reduce lateral melting when ice fraction is low : the purpose is to avoid
385 C disappearing of (up to hThinIce thick) sea-ice by over doing lateral melting
386 C (which would bring iceFrac below iceMaskMin).
387 IF ( iceFrac.LE.lowIcFrac1 ) THEN
388 frace = 0. _d 0
389 ELSEIF (iceFrac.LE.lowIcFrac2 ) THEN
390 frace = MIN( frace, fracEnMelt )
391 ENDIF
392 #endif /* THSICE_FRACEN_POWERLAW */
393
394 c IF (Tsf .EQ. 0. _d 0 .AND. sHeating.GT.0. _d 0) THEN
395 IF ( sHeating.GT.0. _d 0 ) THEN
396 etop = (1. _d 0-frace)*sHeating * dt
397 etope = frace*sHeating * dt
398 ELSE
399 etop = 0. _d 0
400 etope = 0. _d 0
401 C jmc: found few cases where Tsf=0 & sHeating < 0 : add this line to conserv energy:
402 esurp = sHeating * dt
403 ENDIF
404 C-- flux at the base of sea-ice:
405 C conduction H.flx= flxCnB (+ =down); oceanic turbulent H.flx= Fbot (+ =down).
406 C- ==> energy available(+ => melt)= (flxCnB-Fbot)*dt
407 c IF (frzmlt.LT.0. _d 0) THEN
408 c ebot = (1. _d 0-frace)*(flxCnB-Fbot) * dt
409 c ebote = frace*(flxCnB-Fbot) * dt
410 c ELSE
411 c ebot = (flxCnB-Fbot) * dt
412 c ebote = 0. _d 0
413 c ENDIF
414 C- original formulation(above): Loose energy when flxCnB < Fbot < 0
415 ebot = (flxCnB(i,j)-Fbot) * dt
416 IF (ebot.GT.0. _d 0) THEN
417 ebote = frace*ebot
418 ebot = ebot-ebote
419 ELSE
420 ebote = 0. _d 0
421 ENDIF
422 #ifdef ALLOW_DBUG_THSICE
423 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
424 & 'ThSI_CALC_TH: etop,etope,ebot,ebote=', etop,etope,ebot,ebote
425 #endif
426
427 C Initialize layer thicknesses.
428 C Make sure internal ice temperatures do not exceed Tmlt.
429 C If they do, then eliminate the layer. (Dont think this will happen
430 C for reasonable values of i0.)
431
432 hlyr = hi * rec_nlyr
433 DO k = 1, nlyr
434 hnew(k) = hlyr
435 ENDDO
436
437 C Top melt: snow, then ice.
438
439 IF (etop .GT. 0. _d 0) THEN
440 #ifdef ALLOW_AUTODIFF_TAMC
441 CADJ STORE etop = comlev1_thsice_1, key=ikey_1
442 #endif
443 IF (hs. gt. 0. _d 0) THEN
444 rq = rhos * qsnow
445 rqh = rq * hs
446 IF (etop .LT. rqh) THEN
447 hs = hs - etop/rq
448 etop = 0. _d 0
449 ELSE
450 etop = etop - rqh
451 hs = 0. _d 0
452 ENDIF
453 ENDIF
454
455 DO k = 1, nlyr
456 #ifdef ALLOW_AUTODIFF_TAMC
457 ikey_2 = k
458 & + nlyr*(i-1)
459 & + nlyr*sNx*(j-1)
460 & + nlyr*sNx*sNy*act1
461 & + nlyr*sNx*sNy*max1*act2
462 & + nlyr*sNx*sNy*max1*max2*act3
463 & + nlyr*sNx*sNy*max1*max2*max3*act4
464 #endif
465 C--
466 #ifdef ALLOW_AUTODIFF_TAMC
467 CADJ STORE etop = comlev1_thsice_2, key=ikey_2
468 CADJ STORE hnew(k) = comlev1_thsice_2, key=ikey_2
469 #endif
470 IF (etop .GT. 0. _d 0) THEN
471 rq = rhoi * qicen(k)
472 rqh = rq * hnew(k)
473 IF (etop .LT. rqh) THEN
474 hnew(k) = hnew(k) - etop / rq
475 etop = 0. _d 0
476 ELSE
477 etop = etop - rqh
478 hnew(k) = 0. _d 0
479 ENDIF
480 ENDIF
481 ENDDO
482 ELSE
483 etop=0. _d 0
484 ENDIF
485 C If ice is gone and melting energy remains
486 c IF (etop .GT. 0. _d 0) THEN
487 c WRITE (6,*) 'QQ All ice melts from top ', i,j
488 c hi=0. _d 0
489 c go to 200
490 c ENDIF
491
492
493 C Bottom melt/growth.
494
495 IF (ebot .LT. 0. _d 0) THEN
496 C Compute enthalpy of new ice growing at bottom surface.
497 qbot = -cpIce *Tf + Lfresh
498 dhi = -ebot / (qbot * rhoi)
499 ebot = 0. _d 0
500 cph k = nlyr
501 #ifdef ALLOW_AUTODIFF_TAMC
502 CADJ STORE hnew(:) = comlev1_thsice_1, key=ikey_1
503 CADJ STORE qicen(:) = comlev1_thsice_1, key=ikey_1
504 #endif
505 qicen(nlyr) = (hnew(nlyr)*qicen(nlyr)+dhi*qbot) /
506 & (hnew(nlyr)+dhi)
507 hnew(nlyr) = hnew(nlyr) + dhi
508 ELSE
509 DO k = nlyr, 1, -1
510 #ifdef ALLOW_AUTODIFF_TAMC
511 ikey_2 = (nlyr-k+1)
512 & + nlyr*(i-1)
513 & + nlyr*sNx*(j-1)
514 & + nlyr*sNx*sNy*act1
515 & + nlyr*sNx*sNy*max1*act2
516 & + nlyr*sNx*sNy*max1*max2*act3
517 & + nlyr*sNx*sNy*max1*max2*max3*act4
518 #endif
519 C--
520 #ifdef ALLOW_AUTODIFF_TAMC
521 CADJ STORE ebot = comlev1_thsice_2, key=ikey_2
522 CADJ STORE hnew(k) = comlev1_thsice_2, key=ikey_2
523 CADJ STORE qicen(k) = comlev1_thsice_2, key=ikey_2
524 #endif
525 IF (ebot.GT.0. _d 0 .AND. hnew(k).GT.0. _d 0) THEN
526 rq = rhoi * qicen(k)
527 rqh = rq * hnew(k)
528 IF (ebot .LT. rqh) THEN
529 hnew(k) = hnew(k) - ebot / rq
530 ebot = 0. _d 0
531 ELSE
532 ebot = ebot - rqh
533 hnew(k) = 0. _d 0
534 ENDIF
535 ENDIF
536 ENDDO
537
538 C If ice melts completely and snow is left, remove the snow with
539 C energy from the mixed layer
540
541 IF (ebot.GT.0. _d 0 .AND. hs.GT.0. _d 0) THEN
542 rq = rhos * qsnow
543 rqh = rq * hs
544 IF (ebot .LT. rqh) THEN
545 hs = hs - ebot / rq
546 ebot = 0. _d 0
547 ELSE
548 ebot = ebot - rqh
549 hs = 0. _d 0
550 ENDIF
551 ENDIF
552 c IF (ebot .GT. 0. _d 0) THEN
553 c IF (dBug) WRITE(6,*) 'All ice (& snow) melts from bottom'
554 c hi=0. _d 0
555 c go to 200
556 c ENDIF
557 ENDIF
558
559 C Compute new total ice thickness.
560 hi = 0. _d 0
561 DO k = 1, nlyr
562 hi = hi + hnew(k)
563 ENDDO
564 #ifdef ALLOW_DBUG_THSICE
565 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
566 & 'ThSI_CALC_TH: etop, ebot, hi, hs =', etop, ebot, hi, hs
567 #endif
568
569 C If hi < hIceMin, melt the ice.
570 IF ( hi.LT.hIceMin .AND. (hi+hs).GT.0. _d 0 ) THEN
571 esurp = esurp - rhos*qsnow*hs
572 DO k = 1, nlyr
573 esurp = esurp - rhoi*qicen(k)*hnew(k)
574 ENDDO
575 hi = 0. _d 0
576 hs = 0. _d 0
577 Tsf=0. _d 0
578 iceFrac = 0. _d 0
579 DO k = 1, nlyr
580 qicen(k) = 0. _d 0
581 ENDDO
582 #ifdef ALLOW_DBUG_THSICE
583 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
584 & 'ThSI_CALC_TH: -1 : esurp=',esurp
585 #endif
586 ENDIF
587
588 C-- do a mass-budget of sea-ice to compute "fresh" = the fresh-water flux
589 C that is returned to the ocean ; needs to be done before snow/evap
590 fresh = (mwater0 - (rhos*hs + rhoi*hi))/dt
591
592 IF ( hi .LE. 0. _d 0 ) THEN
593 C- return snow to the ocean (account for Latent heat of freezing)
594 fresh = fresh + snowPr
595 qleft = qleft - snowPr*Lfresh
596
597 ELSE
598 C- else: hi > 0
599
600 C Let it snow
601
602 hs = hs + dt*snowPr/rhos
603
604 C If ice evap is used to sublimate surface snow/ice or
605 C if no ice pass on to ocean
606 #ifdef ALLOW_AUTODIFF_TAMC
607 CADJ STORE evap = comlev1_thsice_1, key=ikey_1
608 CADJ STORE hs = comlev1_thsice_1, key=ikey_1
609 #endif
610 IF (hs.GT.0. _d 0) THEN
611 IF (evap/rhos *dt.GT.hs) THEN
612 evap=evap-hs*rhos/dt
613 hs=0. _d 0
614 ELSE
615 hs = hs - evap/rhos *dt
616 evap=0. _d 0
617 ENDIF
618 ENDIF
619 #ifdef ALLOW_AUTODIFF_TAMC
620 CADJ STORE evap = comlev1_thsice_1, key=ikey_1
621 CADJ STORE hi = comlev1_thsice_1, key=ikey_1
622 CADJ STORE hnew(:) = comlev1_thsice_1, key=ikey_1
623 CADJ STORE qicen(:) = comlev1_thsice_1, key=ikey_1
624 #endif
625 IF (hi.GT.0. _d 0.AND.evap.GT.0. _d 0) THEN
626 DO k = 1, nlyr
627 #ifdef ALLOW_AUTODIFF_TAMC
628 ikey_2 = k
629 & + nlyr*(i-1)
630 & + nlyr*sNx*(j-1)
631 & + nlyr*sNx*sNy*act1
632 & + nlyr*sNx*sNy*max1*act2
633 & + nlyr*sNx*sNy*max1*max2*act3
634 & + nlyr*sNx*sNy*max1*max2*max3*act4
635 #endif
636 C--
637 #ifdef ALLOW_AUTODIFF_TAMC
638 CADJ STORE evap = comlev1_thsice_2, key=ikey_2
639 CADJ STORE hnew(k) = comlev1_thsice_2, key=ikey_2
640 CADJ STORE qicen(k) = comlev1_thsice_2, key=ikey_2
641 #endif
642 IF (evap .GT. 0. _d 0) THEN
643 C-- original scheme, does not care about ice temp.
644 C- this can produce small error (< 1.W/m2) in the Energy budget
645 c IF (evap/rhoi *dt.GT.hnew(k)) THEN
646 c evap=evap-hnew(k)*rhoi/dt
647 c hnew(k)=0. _d 0
648 c ELSE
649 c hnew(k) = hnew(k) - evap/rhoi *dt
650 c evap=0. _d 0
651 c ENDIF
652 C-- modified scheme. taking into account Ice enthalpy
653 dhi = evap/rhoi*dt
654 IF (dhi.GE.hnew(k)) THEN
655 evap=evap-hnew(k)*rhoi/dt
656 esurp = esurp - hnew(k)*rhoi*(qicen(k)-Lfresh)
657 hnew(k)=0. _d 0
658 ELSE
659 #ifdef ALLOW_AUTODIFF_TAMC
660 CADJ STORE hnew(k) = comlev1_thsice_2, key=ikey_2
661 #endif
662 hq = hnew(k)*qicen(k)-dhi*Lfresh
663 hnew(k) = hnew(k) - dhi
664 qicen(k)=hq/hnew(k)
665 evap=0. _d 0
666 ENDIF
667 C-------
668 ENDIF
669 ENDDO
670 ENDIF
671 c IF (evap .GT. 0. _d 0) THEN
672 c WRITE (6,*) 'BB All ice sublimates', i,j
673 c hi=0. _d 0
674 c go to 200
675 c ENDIF
676
677 C Compute new total ice thickness.
678
679 hi = 0. _d 0
680 DO k = 1, nlyr
681 hi = hi + hnew(k)
682 ENDDO
683
684 C If hi < hIceMin, melt the ice.
685 IF ( hi.GT.0. _d 0 .AND. hi.LT.hIceMin ) THEN
686 fresh = fresh + (rhos*hs + rhoi*hi)/dt
687 esurp = esurp - rhos*qsnow*hs
688 DO k = 1, nlyr
689 esurp = esurp - rhoi*qicen(k)*hnew(k)
690 ENDDO
691 hi = 0. _d 0
692 hs = 0. _d 0
693 Tsf=0. _d 0
694 iceFrac = 0. _d 0
695 DO k = 1, nlyr
696 qicen(k) = 0. _d 0
697 ENDDO
698 #ifdef ALLOW_DBUG_THSICE
699 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
700 & 'ThSI_CALC_TH: -2 : esurp,fresh=', esurp, fresh
701 #endif
702 ENDIF
703
704 C- else hi > 0: end
705 ENDIF
706
707 IF ( hi .GT. 0. _d 0 ) THEN
708
709 C If there is enough snow to lower the ice/snow interface below
710 C freeboard, convert enough snow to ice to bring the interface back
711 C to sea-level OR if snow height is larger than hsMax, snow is
712 C converted to ice to bring hs down to hsMax. Largest change is
713 C applied and enthalpy of top ice layer adjusted accordingly.
714
715 #ifdef ALLOW_AUTODIFF_TAMC
716 CADJ STORE hi = comlev1_thsice_1, key=ikey_1
717 CADJ STORE hs = comlev1_thsice_1, key=ikey_1
718 CADJ STORE hnew(:) = comlev1_thsice_1, key=ikey_1
719 CADJ STORE qicen(:) = comlev1_thsice_1, key=ikey_1
720 #endif
721 IF ( hs .GT. hi*floodFac .OR. hs .GT. hsMax ) THEN
722 cBB WRITE (6,*) 'Freeboard adjusts'
723 c dhi = (hs * rhos - hi * rhoiw) / rhosw
724 c dhs = dhi * rhoi / rhos
725 dhs = (hs - hi*floodFac) * rhoi / rhosw
726 dhs = MAX( hs - hsMax, dhs )
727 dhi = dhs * rhos / rhoi
728 #ifdef ALLOW_AUTODIFF_TAMC
729 CADJ STORE hnew(:) = comlev1_thsice_1, key=ikey_1
730 #endif
731 rqh = rhoi*qicen(1)*hnew(1) + rhos*qsnow*dhs
732 hnew(1) = hnew(1) + dhi
733 qicen(1) = rqh / (rhoi*hnew(1))
734 hi = hi + dhi
735 hs = hs - dhs
736 ENDIF
737
738
739 C limit ice height
740 C- NOTE: this part does not conserve Energy ;
741 C but surplus of fresh water and salt are taken into account.
742 IF (hi.GT.hiMax) THEN
743 cBB print*,'BBerr, hi>hiMax',i,j,hi
744 chi=hi-hiMax
745 DO k=1,nlyr
746 hnew(k)=hnew(k)-chi/2. _d 0
747 ENDDO
748 fresh = fresh + chi*rhoi/dt
749 ENDIF
750 c IF (hs.GT.hsMax) THEN
751 cc print*,'BBerr, hs>hsMax',i,j,hs
752 c chs=hs-hsMax
753 c hs=hsMax
754 c fresh = fresh + chs*rhos/dt
755 c ENDIF
756
757 C Compute new total ice thickness.
758
759 hi = 0. _d 0
760 DO k = 1, nlyr
761 hi = hi + hnew(k)
762 ENDDO
763
764 #ifdef ALLOW_DBUG_THSICE
765 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
766 & 'ThSI_CALC_TH: b-Winton: hnew, qice =', hnew, qicen
767 #endif
768
769 hlyr = hi * rec_nlyr
770 CALL THSICE_RESHAPE_LAYERS(
771 U qicen,
772 I hlyr, hnew, myThid )
773
774 #ifdef ALLOW_DBUG_THSICE
775 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
776 & 'ThSI_CALC_TH: iceFrac,hi, qtot, hs =', iceFrac,hi,
777 & (qicen(1)+qicen(2))*0.5, hs
778 #endif
779
780 C- if hi > 0 : end
781 ENDIF
782 200 CONTINUE
783
784 C- Compute surplus energy left over from melting.
785
786 IF (hi.LE.0. _d 0) iceFrac=0. _d 0
787
788 C.. heat fluxes left over for ocean
789 qleft = qleft + (Fbot+(esurp+etop+ebot)/dt)
790 #ifdef ALLOW_DBUG_THSICE
791 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
792 & 'ThSI_CALC_TH: [esurp,etop+ebot]/dt =',esurp/dt,etop/dt,ebot/dt
793 #endif
794
795 C-- Evaporation left to the ocean :
796 fresh = fresh - evap
797 C- Correct Atmos. fluxes for this different latent heat:
798 C evap was computed over freezing surf.(Tsf<0), latent heat = Lvap+Lfresh
799 C but should be Lvap only for the fraction "evap" that is left to the ocean.
800 qleft = qleft + evap*Lfresh
801
802 C fresh and salt fluxes
803 c fresh = (mwater0 - (rhos*(hs) + rhoi*(hi)))/dt-evap
804 c fsalt = (msalt0 - rhoi*hi*saltIce)/35. _d 0/dt ! for same units as fresh
805 C note (jmc): fresh is computed from a sea-ice mass budget that already
806 C contains, at this point, snow & evaporation (of snow & ice)
807 C but are not meant to be part of ice/ocean fresh-water flux.
808 C fix: a) like below or b) by making the budget before snow/evap is added
809 c fresh = (mwater0 - (rhos*(hs) + rhoi*(hi)))/dt
810 c & + snow(i,j,bi,bj)*rhos - frwAtm
811 fsalt(i,j) = (msalt0 - rhoi*hi*saltIce)/dt
812
813 #ifdef ALLOW_DBUG_THSICE
814 IF (dBug(i,j,bi,bj) ) THEN
815 WRITE(6,1020)'ThSI_CALC_TH:dH2O,Ev[kg],fresh,fsalt',
816 & (mwater0-(rhos*hs+rhoi*hi))/dt,evap,fresh,fsalt(i,j)
817 WRITE(6,1020)'ThSI_CALC_TH: Qleft,Fbot,extend/dt =',
818 & Qleft,Fbot,(etope+ebote)/dt
819 ENDIF
820 #endif
821
822 C-- add remaining liquid Precip (rain+RunOff) directly to ocean:
823 fresh = fresh + (prcAtm(i,j)-snowPr)
824
825 C-- note: at this point, iceFrac has not been changed (unless reset to zero)
826 C and it can only be reduced by lateral melting in the following part:
827
828 C calculate extent changes
829 extend=etope+ebote
830 IF (iceFrac.GT.0. _d 0.AND.extend.GT.0. _d 0) THEN
831 rq = rhoi * 0.5 _d 0*(qicen(1)+qicen(2))
832 rs = rhos * qsnow
833 rqh = rq * hi + rs * hs
834 freshe=(rhos*hs+rhoi*hi)/dt
835 salte=(rhoi*hi*saltIce)/dt
836 IF ( extend.LT.rqh ) THEN
837 iceFrac=(1. _d 0-extend/rqh)*iceFrac
838 ENDIF
839 IF ( extend.LT.rqh .AND. iceFrac.GE.iceMaskMin ) THEN
840 fresh=fresh+extend/rqh*freshe
841 fsalt(i,j)=fsalt(i,j)+extend/rqh*salte
842 ELSE
843 iceFrac=0. _d 0
844 hi=0. _d 0
845 hs=0. _d 0
846 qleft=qleft+(extend-rqh)/dt
847 fresh=fresh+freshe
848 fsalt(i,j)=fsalt(i,j)+salte
849 ENDIF
850 ELSEIF (extend.GT.0. _d 0) THEN
851 qleft=qleft+extend/dt
852 ENDIF
853
854 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
855 C-- Update output variables :
856
857 C- Diagnostic of Atmos. fresh water flux (E-P) over sea-ice :
858 C substract precip from Evap (<- stored in frwAtm array)
859 frwAtm(i,j) = frwAtm(i,j) - prcAtm(i,j)
860
861 C- update Mixed-Layer Freezing potential heat flux by substracting the
862 C part which has already been accounted for (Fbot):
863 fzMlOc(i,j) = frzmlt - Fbot*iceMask(i,j)
864
865 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
866 #ifdef ALLOW_DBUG_THSICE
867 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
868 & 'ThSI_CALC_TH: iceFrac,flx2oc,fsalt,frw2oc=',
869 & iceFrac, qleft, fsalt(i,j), fresh
870 #endif
871 #ifdef CHECK_ENERGY_CONSERV
872 CALL THSICE_CHECK_CONSERV( dBugFlag, i, j, bi, bj, 0,
873 I iceMask(i,j), iceFrac, hi, hs, qicen,
874 I qleft, fresh, fsalt,
875 I myTime, myIter, myThid )
876 #endif /* CHECK_ENERGY_CONSERV */
877 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
878 C-- Update Sea-Ice state output:
879 icFrac(i,j) = iceFrac
880 hIce(i,j) = hi
881 hSnow(i,j ) = hs
882 tSrf(i,j) = Tsf
883 qIc1(i,j) = qicen(1)
884 qIc2(i,j) = qicen(2)
885 C-- Update oceanic flux output:
886 flx2oc(i,j) = qleft
887 frw2oc(i,j) = fresh
888 c fsalt(i,j) = fsalt
889 ENDIF
890 ENDDO
891 ENDDO
892 #endif /* ALLOW_THSICE */
893
894 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
895
896 RETURN
897 END

  ViewVC Help
Powered by ViewVC 1.1.22