/[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.24 - (show annotations) (download)
Fri Dec 17 04:00:14 2010 UTC (13 years, 4 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62y, checkpoint62x
Changes since 1.23: +23 -23 lines
rename iicekey as ticekey to avoid conflict with pkg/seaice

1 C $Header: /u/gcmpack/MITgcm/pkg/thsice/thsice_calc_thickn.F,v 1.23 2010/10/26 22:26:59 heimbach 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,
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.
30 C.. J. Geophys. 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 "SIZE.h"
67 #include "EEPARAMS.h"
68 #include "THSICE_SIZE.h"
69 #include "THSICE_PARAMS.h"
70 #ifdef ALLOW_AUTODIFF_TAMC
71 # include "tamc.h"
72 # include "tamc_keys.h"
73 #endif
74
75 C !INPUT/OUTPUT PARAMETERS:
76 C == Routine Arguments ==
77 C bi,bj :: tile indices
78 C iMin,iMax :: computation domain: 1rst index range
79 C jMin,jMax :: computation domain: 2nd index range
80 C dBugFlag :: allow to print debugging stuff (e.g. on 1 grid point).
81 C--- Input:
82 C iceMask :: sea-ice fractional mask [0-1]
83 C tFrz :: sea-water freezing temperature [oC] (function of S)
84 C tOce :: surface level oceanic temperature [oC]
85 C v2oc :: square of ocean surface-level velocity [m2/s2]
86 C snowP :: snow precipitation [kg/m2/s]
87 C prcAtm :: total precip from the atmosphere [kg/m2/s]
88 C sHeat :: surf heating flux left to melt snow or ice (= Atmos-conduction)
89 C flxCnB :: heat flux conducted through the ice to bottom surface
90 C--- Modified (input&output):
91 C icFrac :: fraction of grid area covered in ice
92 C hIce :: ice height [m]
93 C hSnow :: snow height [m]
94 C tSrf :: surface (ice or snow) temperature
95 C qIc1 (qicen) :: ice enthalpy (J/kg), 1rst level
96 C qIc2 (qicen) :: ice enthalpy (J/kg), 2nd level
97 C frwAtm (evpAtm):: evaporation to the atmosphere [kg/m2/s] (>0 if evaporate)
98 C fzMlOc :: ocean mixed-layer freezing/melting potential [W/m2]
99 C flx2oc :: net heat flux to ocean [W/m2] (> 0 downward)
100 C--- Output
101 C frw2oc :: Total fresh water flux to ocean [kg/m2/s] (> 0 downward)
102 C fsalt :: salt flux to ocean [g/m2/s] (> 0 downward)
103 C--- Input:
104 C myTime :: current Time of simulation [s]
105 C myIter :: current Iteration number in simulation
106 C myThid :: my Thread Id number
107 INTEGER bi,bj
108 INTEGER iMin, iMax
109 INTEGER jMin, jMax
110 LOGICAL dBugFlag
111 _RL iceMask(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
112 _RL tFrz (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
113 _RL tOce (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
114 _RL v2oc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
115 _RL snowP (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
116 _RL prcAtm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
117 _RL sHeat (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
118 _RL flxCnB (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
119 _RL icFrac (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RL hIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 _RL hSnow (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 _RL tSrf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123 _RL qIc1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124 _RL qIc2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125 _RL frwAtm (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126 _RL fzMlOc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127 _RL flx2oc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128 _RL frw2oc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 _RL fsalt (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 _RL myTime
131 INTEGER myIter
132 INTEGER myThid
133 CEOP
134
135 #ifdef ALLOW_THSICE
136
137 C !LOCAL VARIABLES:
138 C--- local copy of input/output argument list variables (see description above)
139 _RL qicen(1-Olx:sNx+Olx,1-Oly:sNy+Oly,nlyr)
140
141 C == Local Variables ==
142 C i,j,k :: loop indices
143 C rec_nlyr :: reciprocal of number of ice layers (real value)
144 C evapLoc :: evaporation over snow/ice [kg/m2/s] (>0 if evaporate)
145 C Fbot :: oceanic heat flux used to melt/form ice [W/m2]
146 C etop :: energy for top melting (J m-2)
147 C ebot :: energy for bottom melting (J m-2)
148 C etope :: energy (from top) for lateral melting (J m-2)
149 C ebote :: energy (from bottom) for lateral melting (J m-2)
150 C extend :: total energy for lateral melting (J m-2)
151 C hnew(nlyr) :: new ice layer thickness (m)
152 C hlyr :: individual ice layer thickness (m)
153 C dhi :: change in ice thickness
154 C dhs :: change in snow thickness
155 C rq :: rho * q for a layer
156 C rqh :: rho * q * h for a layer
157 C qbot :: enthalpy for new ice at bottom surf (J/kg)
158 C dt :: timestep
159 C esurp :: surplus energy from melting (J m-2)
160 C mwater0 :: fresh water mass gained/lost (kg/m^2)
161 C msalt0 :: salt gained/lost (kg/m^2)
162 C freshe :: fresh water gain from extension melting
163 C salte :: salt gained from extension melting
164 C lowIcFrac1 :: ice-fraction lower limit above which partial (lowIcFrac1)
165 C lowIcFrac2 :: or full (lowIcFrac2) lateral melting is allowed.
166 C from THSICE_RESHAPE_LAYERS
167 C f1 :: Fraction of upper layer ice in new layer
168 C qh1, qh2 :: qice*h for layers 1 and 2
169 C qhtot :: qh1 + qh2
170 C q2tmp :: Temporary value of qice for layer 2
171 INTEGER i,j,k
172 _RL rec_nlyr
173 _RL evapLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
174 _RL Fbot (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
175 _RL etop (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
176 _RL ebot (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
177 _RL etope (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
178 _RL ebote (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
179 _RL esurp (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
180 _RL extend
181 _RL hnew (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nlyr)
182 _RL hlyr
183 _RL dhi
184 _RL dhs
185 _RL rq
186 _RL rqh
187 _RL qbot
188 _RL dt
189 _RL mwater0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
190 _RL msalt0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
191 _RL freshe
192 _RL salte
193 _RL lowIcFrac1, lowIcFrac2
194 _RL f1
195 _RL qh1, qh2
196 _RL qhtot
197 _RL q2tmp
198 #ifdef CHECK_ENERGY_CONSERV
199 _RL qaux(nlyr)
200 #endif /* CHECK_ENERGY_CONSERV */
201
202 _RL ustar, cpchr
203 _RL chi
204 _RL frace, rs, hq
205 #ifdef THSICE_FRACEN_POWERLAW
206 INTEGER powerLaw
207 _RL rec_pLaw
208 _RL c1Mlt, c2Mlt, aMlt, hMlt
209 _RL c1Frz, c2Frz, aFrz, hFrz
210 _RL enFrcMlt(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
211 _RL xxMlt, tmpMlt
212 _RL enFrcFrz(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
213 _RL xxFrz, tmpFrz
214 #endif
215
216 C- define grid-point location where to print debugging values
217 #include "THSICE_DEBUG.h"
218
219 1010 FORMAT(A,I3,3F8.3)
220 1020 FORMAT(A,1P4E11.3)
221
222 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
223
224 #ifdef ALLOW_AUTODIFF_TAMC
225 act1 = bi - myBxLo(myThid)
226 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
227 act2 = bj - myByLo(myThid)
228 max2 = myByHi(myThid) - myByLo(myThid) + 1
229 act3 = myThid - 1
230 max3 = nTx*nTy
231 act4 = ikey_dynamics - 1
232 ticekey = (act1 + 1) + act2*max1
233 & + act3*max1*max2
234 & + act4*max1*max2*max3
235 #endif /* ALLOW_AUTODIFF_TAMC */
236
237 rec_nlyr = nlyr
238 rec_nlyr = 1. _d 0 / rec_nlyr
239 dt = thSIce_deltaT
240
241 C for now, use hard coded threshold (iceMaskMin +1.% and +10.%)
242 lowIcFrac1 = iceMaskMin*1.01 _d 0
243 lowIcFrac2 = iceMaskMin*1.10 _d 0
244 #ifdef THSICE_FRACEN_POWERLAW
245 IF ( powerLawExp2 .GE. 1 ) THEN
246 powerLaw = 1 + 2**powerLawExp2
247 rec_pLaw = powerLaw
248 rec_pLaw = 1. _d 0 / rec_pLaw
249 C- Coef for melting:
250 C lateral-melting energy fraction = fracEnMelt - [ aMlt*(hi-hMlt) ]^powerLaw
251 c1Mlt = fracEnMelt**rec_pLaw
252 c2Mlt = (1. _d 0 - fracEnMelt)**rec_pLaw
253 aMlt = (c1Mlt+c2Mlt)/(hThickIce-hThinIce)
254 hMlt = hThinIce+c2Mlt/aMlt
255 C- Coef for freezing:
256 C thickening energy fraction = fracEnFreez - [ aFrz*(hi-hFrz) ]^powerLaw
257 c1Frz = fracEnFreez**rec_pLaw
258 c2Frz = (1. _d 0 -fracEnFreez)**rec_pLaw
259 aFrz = (c1Frz+c2Frz)/(hThickIce-hThinIce)
260 hFrz = hThinIce+c2Frz/aFrz
261 ELSE
262 C- Linear relation
263 powerLaw = 1
264 aMlt = -1. _d 0 /(hThickIce-hThinIce)
265 hMlt = hThickIce
266 aFrz = -1. _d 0 /(hThickIce-hThinIce)
267 hFrz = hThickIce
268 ENDIF
269 #endif /* THSICE_FRACEN_POWERLAW */
270
271
272 C initialise local arrays
273 DO j=1-Oly,sNy+Oly
274 DO i=1-Olx,sNx+Olx
275 evapLoc(i,j) = 0. _d 0
276 Fbot (i,j) = 0. _d 0
277 etop (i,j) = 0. _d 0
278 ebot (i,j) = 0. _d 0
279 etope (i,j) = 0. _d 0
280 ebote (i,j) = 0. _d 0
281 esurp (i,j) = 0. _d 0
282 mwater0(i,j) = 0. _d 0
283 msalt0 (i,j) = 0. _d 0
284 #ifdef THSICE_FRACEN_POWERLAW
285 enFrcMlt(i,j)= 0. _d 0
286 enFrcFrz(i,j)= 0. _d 0
287 #endif
288 ENDDO
289 ENDDO
290 DO k = 1,nlyr
291 DO j=1-Oly,sNy+Oly
292 DO i=1-Olx,sNx+Olx
293 hnew(i,j,k) = 0. _d 0
294 ENDDO
295 ENDDO
296 ENDDO
297
298 DO j = jMin, jMax
299 DO i = iMin, iMax
300 CML#ifdef ALLOW_AUTODIFF_TAMC
301 CML ikey_1 = i
302 CML & + sNx*(j-1)
303 CML & + sNx*sNy*act1
304 CML & + sNx*sNy*max1*act2
305 CML & + sNx*sNy*max1*max2*act3
306 CML & + sNx*sNy*max1*max2*max3*act4
307 CML#endif /* ALLOW_AUTODIFF_TAMC */
308 CML#ifdef ALLOW_AUTODIFF_TAMC
309 CMLCADJ STORE frwatm(i,j) = comlev1_thsice_1, key=ikey_1
310 CMLCADJ STORE fzmloc(i,j) = comlev1_thsice_1, key=ikey_1
311 CMLCADJ STORE hice(i,j) = comlev1_thsice_1, key=ikey_1
312 CMLCADJ STORE hsnow(i,j) = comlev1_thsice_1, key=ikey_1
313 CMLCADJ STORE icfrac(i,j) = comlev1_thsice_1, key=ikey_1
314 CMLCADJ STORE qic1(i,j) = comlev1_thsice_1, key=ikey_1
315 CMLCADJ STORE qic2(i,j) = comlev1_thsice_1, key=ikey_1
316 CML#endif
317
318 IF (iceMask(i,j).GT.0. _d 0) THEN
319 qicen(i,j,1)= qIc1(i,j)
320 qicen(i,j,2)= qIc2(i,j)
321 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
322 C initialize energies
323 esurp(i,j) = 0. _d 0
324
325 c make a local copy of evaporation
326 evapLoc(i,j) = frwAtm(i,j)
327
328 C------------------------------------------------------------------------
329 C-- Compute growth and/or melting at the top and bottom surfaces
330 C------------------------------------------------------------------------
331
332 #ifdef THSICE_FRACEN_POWERLAW
333 xxMlt = aMlt*(hIce(i,j)-hMlt)
334 xxFrz = aFrz*(hIce(i,j)-hFrz)
335 c--
336 IF ( powerLawExp2 .GE. 1 ) THEN
337 #ifdef TARGET_NEC_SX
338 C avoid the short inner loop that cannot be vectorized
339 xxMlt = xxMlt**powerLaw
340 xxFrz = xxFrz**powerLaw
341 #else
342 tmpMlt = xxMlt
343 tmpFrz = xxFrz
344 DO k=1,powerLawExp2
345 tmpMlt = tmpMlt*tmpMlt
346 tmpFrz = tmpFrz*tmpFrz
347 ENDDO
348 xxMlt = xxMlt*tmpMlt
349 xxFrz = xxFrz*tmpFrz
350 #endif /* TARGET_NEC_SX */
351 xxMlt = fracEnMelt -xxMlt
352 xxFrz = fracEnFreez-xxFrz
353 ENDIF
354 enFrcMlt(i,j) = MAX( 0. _d 0, MIN( xxMlt, 1. _d 0 ) )
355 enFrcFrz(i,j) = MAX( 0. _d 0, MIN( xxFrz, 1. _d 0 ) )
356 #endif /* THSICE_FRACEN_POWERLAW */
357
358 IF (fzMlOc(i,j).GE. 0. _d 0) THEN
359 C !-----------------------------------------------------------------
360 C ! freezing conditions
361 C !-----------------------------------------------------------------
362 Fbot(i,j) = fzMlOc(i,j)
363 IF ( icFrac(i,j).LT.iceMaskMax ) THEN
364 #ifdef THSICE_FRACEN_POWERLAW
365 Fbot(i,j) = enFrcFrz(i,j)*fzMlOc(i,j)
366 #else /* THSICE_FRACEN_POWERLAW */
367 IF (hIce(i,j).GT.hThickIce) THEN
368 C if higher than hThickIce, use all fzMlOc energy to grow extra ice
369 Fbot(i,j) = 0. _d 0
370 ELSEIF (hIce(i,j).GE.hThinIce) THEN
371 C between hThinIce & hThickIce, use partition factor fracEnFreez
372 Fbot(i,j) = (1. _d 0 - fracEnFreez)*fzMlOc(i,j)
373 ENDIF
374 #endif /* THSICE_FRACEN_POWERLAW */
375 ENDIF
376 ELSE
377 C !-----------------------------------------------------------------
378 C ! melting conditions
379 C !-----------------------------------------------------------------
380 C for no currents:
381 ustar = 5. _d -2
382 C frictional velocity between ice and water
383 IF (v2oc(i,j) .NE. 0.)
384 & ustar = SQRT(0.00536 _d 0*v2oc(i,j))
385 ustar=max(5. _d -3,ustar)
386 cpchr =cpWater*rhosw*bMeltCoef
387 Fbot(i,j) = cpchr*(tFrz(i,j)-tOce(i,j))*ustar
388 C fzMlOc < Fbot < 0
389 Fbot(i,j) = max(Fbot(i,j),fzMlOc(i,j))
390 Fbot(i,j) = min(Fbot(i,j),0. _d 0)
391 ENDIF
392
393 C mass of fresh water and salt initially present in ice
394 mwater0(i,j) = rhos*hSnow(i,j) + rhoi*hIce(i,j)
395 msalt0 (i,j) = rhoi*hIce(i,j)*saltIce
396
397 #ifdef ALLOW_DBUG_THSICE
398 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
399 & 'ThSI_CALC_TH: evpAtm, fzMlOc, Fbot =',
400 & frwAtm(i,j),fzMlOc(i,j),Fbot(i,j)
401 #endif
402 C endif iceMask > 0
403 ENDIF
404 C end i/j-loops
405 ENDDO
406 ENDDO
407 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
408 DO j = jMin, jMax
409 DO i = iMin, iMax
410 IF (iceMask(i,j).GT.0. _d 0) THEN
411
412 C Compute energy available for melting/growth.
413
414 #ifdef THSICE_FRACEN_POWERLAW
415 IF ( fracEnMelt.EQ.0. _d 0 ) THEN
416 frace = 0. _d 0
417 ELSE
418 frace = (icFrac(i,j) - lowIcFrac1)/(lowIcFrac2-iceMaskMin)
419 frace = MIN( enFrcMlt(i,j), MAX( 0. _d 0, frace ) )
420 ENDIF
421 #else /* THSICE_FRACEN_POWERLAW */
422 IF ( hIce(i,j).GT.hThickIce .OR. fracEnMelt.EQ.0. _d 0 ) THEN
423 C above certain height (or when no ice fractionation), only melt from top
424 frace = 0. _d 0
425 ELSEIF (hIce(i,j).LT.hThinIce) THEN
426 C below a certain height, all energy goes to changing ice extent
427 frace = 1. _d 0
428 ELSE
429 frace = fracEnMelt
430 ENDIF
431 C Reduce lateral melting when ice fraction is low : the purpose is to avoid
432 C disappearing of (up to hThinIce thick) sea-ice by over doing lateral melting
433 C (which would bring icFrac below iceMaskMin).
434 IF ( icFrac(i,j).LE.lowIcFrac1 ) THEN
435 frace = 0. _d 0
436 ELSEIF (icFrac(i,j).LE.lowIcFrac2 ) THEN
437 frace = MIN( frace, fracEnMelt )
438 ENDIF
439 #endif /* THSICE_FRACEN_POWERLAW */
440
441 c IF (tSrf(i,j) .EQ. 0. _d 0 .AND. sHeat(i,j).GT.0. _d 0) THEN
442 IF ( sHeat(i,j).GT.0. _d 0 ) THEN
443 etop(i,j) = (1. _d 0-frace)*sHeat(i,j) * dt
444 etope(i,j) = frace*sHeat(i,j) * dt
445 ELSE
446 etop(i,j) = 0. _d 0
447 etope(i,j) = 0. _d 0
448 C jmc: found few cases where tSrf=0 & sHeat < 0 : add this line to conserv energy:
449 esurp(i,j) = sHeat(i,j) * dt
450 ENDIF
451 C-- flux at the base of sea-ice:
452 C conduction H.flx= flxCnB (+ =down); oceanic turbulent H.flx= Fbot (+ =down).
453 C- ==> energy available(+ => melt)= (flxCnB-Fbot)*dt
454 c IF (fzMlOc(i,j).LT.0. _d 0) THEN
455 c ebot(i,j) = (1. _d 0-frace)*(flxCnB-Fbot(i,j)) * dt
456 c ebote(i,j) = frace*(flxCnB-Fbot(i,j)) * dt
457 c ELSE
458 c ebot(i,j) = (flxCnB-Fbot(i,j)) * dt
459 c ebote(i,j) = 0. _d 0
460 c ENDIF
461 C- original formulation(above): Loose energy when flxCnB < Fbot < 0
462 ebot(i,j) = (flxCnB(i,j)-Fbot(i,j)) * dt
463 IF (ebot(i,j).GT.0. _d 0) THEN
464 ebote(i,j) = frace*ebot(i,j)
465 ebot(i,j) = ebot(i,j)-ebote(i,j)
466 ELSE
467 ebote(i,j) = 0. _d 0
468 ENDIF
469 #ifdef ALLOW_DBUG_THSICE
470 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
471 & 'ThSI_CALC_TH: etop,etope,ebot,ebote=',
472 & etop(i,j),etope(i,j),ebot(i,j),ebote(i,j)
473 #endif
474 C endif iceMask > 0
475 ENDIF
476 C end i/j-loops
477 ENDDO
478 ENDDO
479
480 C Initialize layer thicknesses. Divide total thickness equally between
481 C layers
482 DO k = 1, nlyr
483 DO j = jMin, jMax
484 DO i = iMin, iMax
485 hnew(i,j,k) = hIce(i,j) * rec_nlyr
486 ENDDO
487 ENDDO
488 ENDDO
489
490 DO j = jMin, jMax
491 DO i = iMin, iMax
492 CML#ifdef ALLOW_AUTODIFF_TAMC
493 CMLCADJ STORE etop(i,j) = comlev1_thsice_1, key=ikey_1
494 CML#endif
495 IF (iceMask(i,j) .GT. 0. _d 0 .AND.
496 & etop(i,j) .GT. 0. _d 0 .AND.
497 & hSnow(i,j) .GT. 0. _d 0) THEN
498
499 C Make sure internal ice temperatures do not exceed Tmlt.
500 C If they do, then eliminate the layer. (Dont think this will happen
501 C for reasonable values of i0.)
502 C Top melt: snow, then ice.
503 rq = rhos * qsnow
504 rqh = rq * hSnow(i,j)
505 IF (etop(i,j) .LT. rqh) THEN
506 hSnow(i,j) = hSnow(i,j) - etop(i,j)/rq
507 etop(i,j) = 0. _d 0
508 ELSE
509 etop(i,j) = etop(i,j) - rqh
510 hSnow(i,j) = 0. _d 0
511 ENDIF
512 C endif iceMask > 0, etc.
513 ENDIF
514 C end i/j-loops
515 ENDDO
516 ENDDO
517 C two layers of ice
518 DO k = 1, nlyr
519 DO j = jMin, jMax
520 DO i = iMin, iMax
521 IF (iceMask(i,j).GT.0. _d 0) THEN
522 CML#ifdef ALLOW_AUTODIFF_TAMC
523 CML ikey_2 = k
524 CML & + nlyr*(i-1)
525 CML & + nlyr*sNx*(j-1)
526 CML & + nlyr*sNx*sNy*act1
527 CML & + nlyr*sNx*sNy*max1*act2
528 CML & + nlyr*sNx*sNy*max1*max2*act3
529 CML & + nlyr*sNx*sNy*max1*max2*max3*act4
530 CML#endif
531 CML#ifdef ALLOW_AUTODIFF_TAMC
532 CMLCADJ STORE etop(i,j) = comlev1_thsice_2, key=ikey_2
533 CMLCADJ STORE hnew(i,j,k) = comlev1_thsice_2, key=ikey_2
534 CML#endif
535 IF (etop(i,j) .GT. 0. _d 0) THEN
536 rq = rhoi * qicen(i,j,k)
537 rqh = rq * hnew(i,j,k)
538 IF (etop(i,j) .LT. rqh) THEN
539 hnew(i,j,k) = hnew(i,j,k) - etop(i,j) / rq
540 etop(i,j) = 0. _d 0
541 ELSE
542 etop(i,j) = etop(i,j) - rqh
543 hnew(i,j,k) = 0. _d 0
544 ENDIF
545 ELSE
546 etop(i,j)=0. _d 0
547 ENDIF
548 C If ice is gone and melting energy remains
549 c IF (etop(i,j) .GT. 0. _d 0) THEN
550 c WRITE (6,*) 'QQ All ice melts from top ', i,j
551 c hIce(i,j)=0. _d 0
552 c go to 200
553 c ENDIF
554
555 C endif iceMask > 0
556 ENDIF
557 C end i/j-loops
558 ENDDO
559 ENDDO
560 C end k-loop
561 ENDDO
562
563 DO j = jMin, jMax
564 DO i = iMin, iMax
565 IF (iceMask(i,j).GT.0. _d 0 .AND. ebot(i,j) .LT. 0. _d 0) THEN
566 C Bottom melt/growth.
567 C Compute enthalpy of new ice growing at bottom surface.
568 qbot = -cpIce *tFrz(i,j) + Lfresh
569 dhi = -ebot(i,j) / (qbot * rhoi)
570 ebot(i,j) = 0. _d 0
571 cph k = nlyr
572 CML#ifdef ALLOW_AUTODIFF_TAMC
573 CMLCADJ STORE hnew(i,j,:) = comlev1_thsice_1, key=ikey_1
574 CMLCADJ STORE qicen(i,j,:) = comlev1_thsice_1, key=ikey_1
575 CML#endif
576 qicen(i,j,nlyr) =
577 & (hnew(i,j,nlyr)*qicen(i,j,nlyr)+dhi*qbot) /
578 & (hnew(i,j,nlyr)+dhi)
579 hnew(i,j,nlyr) = hnew(i,j,nlyr) + dhi
580
581 C endif iceMask > 0 and ebot < 0
582 ENDIF
583 C end i/j-loops
584 ENDDO
585 ENDDO
586
587 #ifdef ALLOW_AUTODIFF_TAMC
588 CADJ STORE etop(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
589 CADJ STORE ebot(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
590 CADJ STORE hnew(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
591 CADJ STORE qicen(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
592 #endif
593
594 DO k = nlyr, 1, -1
595 CML#ifdef ALLOW_AUTODIFF_TAMC
596 CML ikey_2 = (nlyr-k+1)
597 CML & + nlyr*(i-1)
598 CML & + nlyr*sNx*(j-1)
599 CML & + nlyr*sNx*sNy*act1
600 CML & + nlyr*sNx*sNy*max1*act2
601 CML & + nlyr*sNx*sNy*max1*max2*act3
602 CML & + nlyr*sNx*sNy*max1*max2*max3*act4
603 CML#endif
604 CML#ifdef ALLOW_AUTODIFF_TAMC
605 CMLCADJ STORE ebot(i,j) = comlev1_thsice_2, key=ikey_2
606 CMLCADJ STORE hnew(i,j,k) = comlev1_thsice_2, key=ikey_2
607 CMLCADJ STORE qicen(i,j,k) = comlev1_thsice_2, key=ikey_2
608 CML#endif
609 DO j = jMin, jMax
610 DO i = iMin, iMax
611 IF (iceMask(i,j) .GT. 0. _d 0 .AND.
612 & ebot(i,j) .GT. 0. _d 0 .AND.
613 & hnew(i,j,k) .GT. 0. _d 0) THEN
614 rq = rhoi * qicen(i,j,k)
615 rqh = rq * hnew(i,j,k)
616 IF (ebot(i,j) .LT. rqh) THEN
617 hnew(i,j,k) = hnew(i,j,k) - ebot(i,j) / rq
618 ebot(i,j) = 0. _d 0
619 ELSE
620 ebot(i,j) = ebot(i,j) - rqh
621 hnew(i,j,k) = 0. _d 0
622 ENDIF
623 C endif iceMask > 0 etc.
624 ENDIF
625 C end i/j-loops
626 ENDDO
627 ENDDO
628 C end k-loop
629 ENDDO
630 C If ice melts completely and snow is left, remove the snow with
631 C energy from the mixed layer
632 DO j = jMin, jMax
633 DO i = iMin, iMax
634 IF (iceMask(i,j) .GT. 0. _d 0 .AND.
635 & ebot(i,j) .GT. 0. _d 0 .AND.
636 & hSnow(i,j) .GT. 0. _d 0) THEN
637 rq = rhos * qsnow
638 rqh = rq * hSnow(i,j)
639 IF (ebot(i,j) .LT. rqh) THEN
640 hSnow(i,j) = hSnow(i,j) - ebot(i,j) / rq
641 ebot(i,j) = 0. _d 0
642 ELSE
643 ebot(i,j) = ebot(i,j) - rqh
644 hSnow(i,j) = 0. _d 0
645 ENDIF
646 c IF (ebot(i,j) .GT. 0. _d 0) THEN
647 c IF (dBug) WRITE(6,*) 'All ice (& snow) melts from bottom'
648 c hIce(i,j)=0. _d 0
649 c go to 200
650 c ENDIF
651
652 C endif iceMask > 0, etc.
653 ENDIF
654 C end i/j-loops
655 ENDDO
656 ENDDO
657 DO j = jMin, jMax
658 DO i = iMin, iMax
659 IF (iceMask(i,j).GT.0. _d 0) THEN
660 C Compute new total ice thickness.
661 hIce(i,j) = hnew(i,j,1) + hnew(i,j,2)
662 #ifdef ALLOW_DBUG_THSICE
663 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
664 & 'ThSI_CALC_TH: etop, ebot, hIce, hSnow =',
665 & etop(i,j), ebot(i,j), hIce(i,j), hSnow(i,j)
666 #endif
667
668 C If hIce < hIceMin, melt the ice.
669 IF ( hIce(i,j).LT.hIceMin
670 & .AND. (hIce(i,j)+hSnow(i,j)).GT.0. _d 0 ) THEN
671 esurp(i,j) = esurp(i,j) - rhos*qsnow*hSnow(i,j)
672 & - rhoi*qicen(i,j,1)*hnew(i,j,1)
673 & - rhoi*qicen(i,j,2)*hnew(i,j,2)
674 hIce(i,j) = 0. _d 0
675 hSnow(i,j) = 0. _d 0
676 tSrf(i,j) = 0. _d 0
677 icFrac(i,j) = 0. _d 0
678 qicen(i,j,1) = 0. _d 0
679 qicen(i,j,2) = 0. _d 0
680 #ifdef ALLOW_DBUG_THSICE
681 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
682 & 'ThSI_CALC_TH: -1 : esurp=',esurp(i,j)
683 #endif
684 ENDIF
685
686 C endif iceMask > 0
687 ENDIF
688 C end i/j-loops
689 ENDDO
690 ENDDO
691
692 DO j = jMin, jMax
693 DO i = iMin, iMax
694 IF (iceMask(i,j).GT.0. _d 0) THEN
695
696 C-- do a mass-budget of sea-ice to compute "fresh" = the fresh-water flux
697 C that is returned to the ocean ; needs to be done before snow/evap
698 frw2oc(i,j) = (mwater0(i,j)
699 & - (rhos*hSnow(i,j)+rhoi*hIce(i,j)))/dt
700
701 IF ( hIce(i,j) .LE. 0. _d 0 ) THEN
702 C- return snow to the ocean (account for Latent heat of freezing)
703 frw2oc(i,j) = frw2oc(i,j) + snowP(i,j)
704 flx2oc(i,j) = flx2oc(i,j) - snowP(i,j)*Lfresh
705 ENDIF
706
707 C endif iceMask > 0
708 ENDIF
709 C end i/j-loops
710 ENDDO
711 ENDDO
712 C- else: hIce > 0
713 DO j = jMin, jMax
714 DO i = iMin, iMax
715 IF (iceMask(i,j).GT.0. _d 0) THEN
716
717 IF ( hIce(i,j) .GT. 0. _d 0 ) THEN
718 C Let it snow
719 hSnow(i,j) = hSnow(i,j) + dt*snowP(i,j)/rhos
720 C If ice evap is used to sublimate surface snow/ice or
721 C if no ice pass on to ocean
722 CML#ifdef ALLOW_AUTODIFF_TAMC
723 CMLCADJ STORE evapLoc(i,j) = comlev1_thsice_1, key=ikey_1
724 CMLCADJ STORE hSnow(i,j) = comlev1_thsice_1, key=ikey_1
725 CML#endif
726 IF (hSnow(i,j).GT.0. _d 0) THEN
727 IF (evapLoc(i,j)/rhos *dt.GT.hSnow(i,j)) THEN
728 evapLoc(i,j)=evapLoc(i,j)-hSnow(i,j)*rhos/dt
729 hSnow(i,j)=0. _d 0
730 ELSE
731 hSnow(i,j) = hSnow(i,j) - evapLoc(i,j)/rhos *dt
732 evapLoc(i,j)=0. _d 0
733 ENDIF
734 ENDIF
735 C endif hice > 0
736 ENDIF
737 C endif iceMask > 0
738 ENDIF
739 C end i/j-loops
740 ENDDO
741 ENDDO
742
743 #ifdef ALLOW_AUTODIFF_TAMC
744 CADJ STORE evaploc(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
745 CADJ STORE hnew(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
746 #endif
747
748 C- else: hIce > 0
749 DO k = 1, nlyr
750 DO j = jMin, jMax
751 DO i = iMin, iMax
752 IF (iceMask(i,j).GT.0. _d 0 ) THEN
753
754 CML#ifdef ALLOW_AUTODIFF_TAMC
755 CMLCADJ STORE evapLoc(i,j) = comlev1_thsice_1, key=ikey_1
756 CMLCADJ STORE hIce(i,j) = comlev1_thsice_1, key=ikey_1
757 CMLCADJ STORE hnew(i,j,:) = comlev1_thsice_1, key=ikey_1
758 CMLCADJ STORE qicen(i,j,:) = comlev1_thsice_1, key=ikey_1
759 CML#endif
760 IF (hIce(i,j).GT.0. _d 0.AND.evapLoc(i,j).GT.0. _d 0) THEN
761 CML#ifdef ALLOW_AUTODIFF_TAMC
762 CML ikey_2 = k
763 CML & + nlyr*(i-1)
764 CML & + nlyr*sNx*(j-1)
765 CML & + nlyr*sNx*sNy*act1
766 CML & + nlyr*sNx*sNy*max1*act2
767 CML & + nlyr*sNx*sNy*max1*max2*act3
768 CML & + nlyr*sNx*sNy*max1*max2*max3*act4
769 CML#endif
770 CMLC--
771 CML#ifdef ALLOW_AUTODIFF_TAMC
772 CMLCADJ STORE evapLoc(i,j) = comlev1_thsice_2, key=ikey_2
773 CMLCADJ STORE hnew(i,j,k) = comlev1_thsice_2, key=ikey_2
774 CMLCADJ STORE qicen(i,j,k) = comlev1_thsice_2, key=ikey_2
775 CML#endif
776 C IF (evapLoc(i,j) .GT. 0. _d 0) THEN
777 C-- original scheme, does not care about ice temp.
778 C- this can produce small error (< 1.W/m2) in the Energy budget
779 c IF (evapLoc(i,j)/rhoi *dt.GT.hnew(i,j,k)) THEN
780 c evapLoc(i,j)=evapLoc(i,j)-hnew(i,j,k)*rhoi/dt
781 c hnew(i,j,k)=0. _d 0
782 c ELSE
783 c hnew(i,j,k) = hnew(i,j,k) - evapLoc(i,j)/rhoi *dt
784 c evapLoc(i,j)=0. _d 0
785 c ENDIF
786 C-- modified scheme. taking into account Ice enthalpy
787 dhi = evapLoc(i,j)/rhoi*dt
788 IF (dhi.GE.hnew(i,j,k)) THEN
789 evapLoc(i,j)=evapLoc(i,j)-hnew(i,j,k)*rhoi/dt
790 esurp(i,j) = esurp(i,j)
791 & - hnew(i,j,k)*rhoi*(qicen(i,j,k)-Lfresh)
792 hnew(i,j,k)=0. _d 0
793 ELSE
794 CML#ifdef ALLOW_AUTODIFF_TAMC
795 CMLCADJ STORE hnew(i,j,k) = comlev1_thsice_2, key=ikey_2
796 CML#endif
797 hq = hnew(i,j,k)*qicen(i,j,k)-dhi*Lfresh
798 hnew(i,j,k) = hnew(i,j,1) - dhi
799 qicen(i,j,k)=hq/hnew(i,j,k)
800 evapLoc(i,j)=0. _d 0
801 ENDIF
802 C-------
803 c IF (evapLoc(i,j) .GT. 0. _d 0) THEN
804 c WRITE (6,*) 'BB All ice sublimates', i,j
805 c hIce(i,j)=0. _d 0
806 c go to 200
807 c ENDIF
808 C endif hice > 0 and evaploc > 0
809 ENDIF
810 C endif iceMask > 0
811 ENDIF
812 C end i/j-loops
813 ENDDO
814 ENDDO
815 C end k-loop
816 ENDDO
817
818 #ifdef ALLOW_AUTODIFF_TAMC
819 CADJ STORE etop(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
820 CADJ STORE icemask(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
821 CADJ STORE hice(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
822 CADJ STORE hnew(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
823 CADJ STORE qicen(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
824 #endif
825
826 C still else: hice > 0
827 DO j = jMin, jMax
828 DO i = iMin, iMax
829 IF (iceMask(i,j).GT.0. _d 0) THEN
830 IF (hIce(i,j) .GT. 0. _d 0) THEN
831 C Compute new total ice thickness.
832 hIce(i,j) = hnew(i,j,1) + hnew(i,j,2)
833 C If hIce < hIceMin, melt the ice.
834 IF ( hIce(i,j).GT.0. _d 0 .AND. hIce(i,j).LT.hIceMin ) THEN
835 frw2oc(i,j) = frw2oc(i,j)
836 & + (rhos*hSnow(i,j) + rhoi*hIce(i,j))/dt
837 esurp(i,j) = esurp(i,j) - rhos*qsnow*hSnow(i,j)
838 & - rhoi*qicen(i,j,1)*hnew(i,j,1)
839 & - rhoi*qicen(i,j,2)*hnew(i,j,2)
840 hIce(i,j) = 0. _d 0
841 hSnow(i,j) = 0. _d 0
842 tSrf(i,j) = 0. _d 0
843 icFrac(i,j) = 0. _d 0
844 qicen(i,j,1) = 0. _d 0
845 qicen(i,j,2) = 0. _d 0
846 #ifdef ALLOW_DBUG_THSICE
847 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
848 & 'ThSI_CALC_TH: -2 : esurp,frw2oc=',
849 & esurp(i,j), frw2oc(i,j)
850 #endif
851 ENDIF
852
853 C-- else hIce > 0: end
854 ENDIF
855
856 C endif iceMask > 0
857 ENDIF
858 C end i/j-loops
859 ENDDO
860 ENDDO
861
862 #ifdef ALLOW_AUTODIFF_TAMC
863 CADJ STORE icemask(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
864 CADJ STORE hice(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
865 CADJ STORE hnew(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
866 CADJ STORE hsnow(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
867 CADJ STORE qicen(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
868 #endif
869
870 DO j = jMin, jMax
871 DO i = iMin, iMax
872 IF (iceMask(i,j).GT.0. _d 0) THEN
873
874 IF ( hIce(i,j) .GT. 0. _d 0 ) THEN
875
876 C If there is enough snow to lower the ice/snow interface below
877 C freeboard, convert enough snow to ice to bring the interface back
878 C to sea-level OR if snow height is larger than hsMax, snow is
879 C converted to ice to bring hSnow down to hsMax. Largest change is
880 C applied and enthalpy of top ice layer adjusted accordingly.
881
882 #ifdef ALLOW_AUTODIFF_TAMC
883 ikey_1 = i
884 & + sNx*(j-1)
885 & + sNx*sNy*act1
886 & + sNx*sNy*max1*act2
887 & + sNx*sNy*max1*max2*act3
888 & + sNx*sNy*max1*max2*max3*act4
889 #endif /* ALLOW_AUTODIFF_TAMC */
890
891 #ifdef ALLOW_AUTODIFF_TAMC
892 CADJ STORE hIce(i,j) = comlev1_thsice_1, key=ikey_1
893 CADJ STORE hSnow(i,j) = comlev1_thsice_1, key=ikey_1
894 CADJ STORE hnew(i,j,:) = comlev1_thsice_1, key=ikey_1
895 CADJ STORE qicen(i,j,:) = comlev1_thsice_1, key=ikey_1
896 #endif
897 IF ( hSnow(i,j) .GT. hIce(i,j)*floodFac
898 & .OR. hSnow(i,j) .GT. hsMax ) THEN
899 cBB WRITE (6,*) 'Freeboard adjusts'
900 c dhi = (hSnow(i,j) * rhos - hIce(i,j) * rhoiw) / rhosw
901 c dhs = dhi * rhoi / rhos
902 #ifdef ALLOW_AUTODIFF_TAMC
903 CADJ STORE hnew(i,j,:) = comlev1_thsice_1, key=ikey_1
904 #endif
905 dhs = (hSnow(i,j) - hIce(i,j)*floodFac) * rhoi / rhosw
906 dhs = MAX( hSnow(i,j) - hsMax, dhs )
907 dhi = dhs * rhos / rhoi
908 rqh = rhoi*qicen(i,j,1)*hnew(i,j,1) + rhos*qsnow*dhs
909 hnew(i,j,1) = hnew(i,j,1) + dhi
910 qicen(i,j,1) = rqh / (rhoi*hnew(i,j,1))
911 hIce(i,j) = hIce(i,j) + dhi
912 hSnow(i,j) = hSnow(i,j) - dhs
913 ENDIF
914
915
916 C limit ice height
917 C- NOTE: this part does not conserve Energy ;
918 C but surplus of fresh water and salt are taken into account.
919 IF (hIce(i,j).GT.hiMax) THEN
920 cBB print*,'BBerr, hIce>hiMax',i,j,hIce(i,j)
921 chi=hIce(i,j)-hiMax
922 hnew(i,j,1)=hnew(i,j,1)-chi/2. _d 0
923 hnew(i,j,2)=hnew(i,j,2)-chi/2. _d 0
924 frw2oc(i,j) = frw2oc(i,j) + chi*rhoi/dt
925 ENDIF
926 c IF (hSnow(i,j).GT.hsMax) THEN
927 cc print*,'BBerr, hSnow>hsMax',i,j,hSnow(i,j)
928 c chs=hSnow(i,j)-hsMax
929 c hSnow(i,j)=hsMax
930 c frw2oc(i,j) = frw2oc(i,j) + chs*rhos/dt
931 c ENDIF
932
933 C Compute new total ice thickness.
934 hIce(i,j) = hnew(i,j,1) + hnew(i,j,2)
935
936 #ifdef ALLOW_DBUG_THSICE
937 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
938 & 'ThSI_CALC_TH: b-Winton: hnew, qice =',
939 & hnew(i,j,1), hnew(i,j,2),
940 & qicen(i,j,1), qicen(i,j,2)
941 #endif
942
943 hlyr = hIce(i,j) * rec_nlyr
944 CML CALL THSICE_RESHAPE_LAYERS(
945 CML U qicen(i,j,:),
946 CML I hlyr, hnew(i,j,:), myThid )
947 C inlined version of S/R THSICE_RESHAPE_LAYERS
948 C | Repartition into equal-thickness layers, conserving energy.
949 C *==========================================================*
950 C | This is the 2-layer version (formerly "NEW_LAYERS_WINTON")
951 C | from M. Winton 1999, JAOT, sea-ice model.
952 if (hnew(i,j,1).gt.hnew(i,j,2)) then
953 C-- Layer 1 gives ice to layer 2
954 f1 = (hnew(i,j,1)-hlyr)/hlyr
955 q2tmp = f1*qicen(i,j,1) + (1. _d 0-f1)*qicen(i,j,2)
956 if (q2tmp.gt.Lfresh) then
957 qicen(i,j,2) = q2tmp
958 else
959 C- Keep q2 fixed to avoid q2<Lfresh and T2>0
960 qh2 = hlyr*qicen(i,j,2)
961 qhtot = hnew(i,j,1)*qicen(i,j,1) + hnew(i,j,2)*qicen(i,j,2)
962 qh1 = qhtot - qh2
963 qicen(i,j,1) = qh1/hlyr
964 endif
965 else
966 C- Layer 2 gives ice to layer 1
967 f1 = hnew(i,j,1)/hlyr
968 qicen(i,j,1) = f1*qicen(i,j,1) + (1. _d 0-f1)*qicen(i,j,2)
969 endif
970 C end of inlined S/R THSICE_RESHAPE_LAYERS
971
972 #ifdef ALLOW_DBUG_THSICE
973 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
974 & 'ThSI_CALC_TH: icFrac,hIce, qtot, hSnow =',
975 & icFrac(i,j),hIce(i,j), (qicen(i,j,1)+qicen(i,j,2))*0.5,
976 & hSnow(i,j)
977 #endif
978
979 C- if hIce > 0 : end
980 ENDIF
981 200 CONTINUE
982
983 C endif iceMask > 0
984 ENDIF
985 C end i/j-loops
986 ENDDO
987 ENDDO
988
989 DO j = jMin, jMax
990 DO i = iMin, iMax
991 IF (iceMask(i,j).GT.0. _d 0) THEN
992 C- Compute surplus energy left over from melting.
993
994 IF (hIce(i,j).LE.0. _d 0) icFrac(i,j)=0. _d 0
995
996 C.. heat fluxes left over for ocean
997 flx2oc(i,j) = flx2oc(i,j)
998 & + (Fbot(i,j)+(esurp(i,j)+etop(i,j)+ebot(i,j))/dt)
999 #ifdef ALLOW_DBUG_THSICE
1000 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
1001 & 'ThSI_CALC_TH: [esurp,etop+ebot]/dt =',
1002 & esurp(i,j)/dt,etop(i,j)/dt,ebot(i,j)/dt
1003 #endif
1004
1005 C-- Evaporation left to the ocean :
1006 frw2oc(i,j) = frw2oc(i,j) - evapLoc(i,j)
1007 C-- Correct Atmos. fluxes for this different latent heat:
1008 C evap was computed over freezing surf.(tSrf<0), latent heat = Lvap+Lfresh
1009 C but should be Lvap only for the fraction "evap" that is left to the ocean.
1010 flx2oc(i,j) = flx2oc(i,j) + evapLoc(i,j)*Lfresh
1011
1012 C fresh and salt fluxes
1013 c frw2oc(i,j) = (mwater0(i,j) - (rhos*(hSnow(i,j))
1014 c & + rhoi*(hIce(i,j))))/dt-evapLoc(i,j)
1015 c fsalt = (msalt0(i,j) - rhoi*hIce(i,j)*saltIce)/35. _d 0/dt ! for same units as frw2oc
1016 C note (jmc): frw2oc is computed from a sea-ice mass budget that already
1017 C contains, at this point, snow & evaporation (of snow & ice)
1018 C but are not meant to be part of ice/ocean fresh-water flux.
1019 C fix: a) like below or b) by making the budget before snow/evap is added
1020 c frw2oc(i,j) = (mwater0(i,j) - (rhos*(hSnow(i,j)) + rhoi*(hIce(i,j))))/dt
1021 c & + snow(i,j,bi,bj)*rhos - frwAtm
1022 fsalt(i,j) = (msalt0(i,j) - rhoi*hIce(i,j)*saltIce)/dt
1023
1024 #ifdef ALLOW_DBUG_THSICE
1025 IF (dBug(i,j,bi,bj) ) THEN
1026 WRITE(6,1020)'ThSI_CALC_TH:dH2O,Ev[kg],frw2oc,fsalt',
1027 & (mwater0(i,j)-(rhos*hSnow(i,j)+rhoi*hIce(i,j)))/dt,
1028 & evapLoc(i,j),frw2oc(i,j),fsalt(i,j)
1029 WRITE(6,1020)'ThSI_CALC_TH: flx2oc,Fbot,extend/dt =',
1030 & flx2oc(I,J),Fbot(i,j),(etope(i,j)+ebote(i,j))/dt
1031 ENDIF
1032 #endif
1033
1034 C-- add remaining liquid Precip (rain+RunOff) directly to ocean:
1035 frw2oc(i,j) = frw2oc(i,j) + (prcAtm(i,j)-snowP(i,j))
1036
1037 C endif iceMask > 0
1038 ENDIF
1039 C end i/j-loops
1040 ENDDO
1041 ENDDO
1042
1043 #ifdef ALLOW_AUTODIFF_TAMC
1044 CADJ STORE icemask(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
1045 CADJ STORE icfrac(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
1046 CADJ STORE hnew(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
1047 CADJ STORE hsnow(:,:) = comlev1_bibj, key=ticekey, byte=isbyte
1048 CADJ STORE qicen(:,:,:) = comlev1_bibj, key=ticekey, byte=isbyte
1049 #endif
1050
1051 DO j = jMin, jMax
1052 DO i = iMin, iMax
1053 IF (iceMask(i,j).GT.0. _d 0) THEN
1054 C-- note: at this point, icFrac has not been changed (unless reset to zero)
1055 C and it can only be reduced by lateral melting in the following part:
1056
1057 C calculate extent changes
1058 extend=etope(i,j)+ebote(i,j)
1059 IF (icFrac(i,j).GT.0. _d 0.AND.extend.GT.0. _d 0) THEN
1060 rq = rhoi * 0.5 _d 0*(qicen(i,j,1)+qicen(i,j,2))
1061 rs = rhos * qsnow
1062 rqh = rq * hIce(i,j) + rs * hSnow(i,j)
1063 freshe=(rhos*hSnow(i,j)+rhoi*hIce(i,j))/dt
1064 salte=(rhoi*hIce(i,j)*saltIce)/dt
1065 IF ( extend.LT.rqh ) THEN
1066 icFrac(i,j)=(1. _d 0-extend/rqh)*icFrac(i,j)
1067 ENDIF
1068 IF ( extend.LT.rqh .AND. icFrac(i,j).GE.iceMaskMin ) THEN
1069 frw2oc(i,j)=frw2oc(i,j)+extend/rqh*freshe
1070 fsalt(i,j)=fsalt(i,j)+extend/rqh*salte
1071 ELSE
1072 icFrac(i,j)=0. _d 0
1073 hIce(i,j) =0. _d 0
1074 hSnow(i,j) =0. _d 0
1075 flx2oc(i,j)=flx2oc(i,j)+(extend-rqh)/dt
1076 frw2oc(i,j)=frw2oc(i,j)+freshe
1077 fsalt(i,j)=fsalt(i,j)+salte
1078 ENDIF
1079 ELSEIF (extend.GT.0. _d 0) THEN
1080 flx2oc(i,j)=flx2oc(i,j)+extend/dt
1081 ENDIF
1082 C endif iceMask > 0
1083 ENDIF
1084 C end i/j-loops
1085 ENDDO
1086 ENDDO
1087 DO j = jMin, jMax
1088 DO i = iMin, iMax
1089 IF (iceMask(i,j).GT.0. _d 0) THEN
1090 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1091 C-- Update output variables :
1092
1093 C-- Diagnostic of Atmos. fresh water flux (E-P) over sea-ice :
1094 C substract precip from Evap (<- stored in frwAtm array)
1095 frwAtm(i,j) = frwAtm(i,j) - prcAtm(i,j)
1096
1097 C-- update Mixed-Layer Freezing potential heat flux by substracting the
1098 C part which has already been accounted for (Fbot):
1099 fzMlOc(i,j) = fzMlOc(i,j) - Fbot(i,j)*iceMask(i,j)
1100
1101 C-- Update Sea-Ice state output:
1102 qIc1(i,j) = qicen(i,j,1)
1103 qIc2(i,j) = qicen(i,j,2)
1104 #ifdef ALLOW_DBUG_THSICE
1105 IF (dBug(i,j,bi,bj) ) WRITE(6,1020)
1106 & 'ThSI_CALC_TH: icFrac,flx2oc,fsalt,frw2oc=',
1107 & icFrac(i,j), flx2oc(i,j), fsalt(i,j), frw2oc(i,j)
1108 #endif
1109 C endif iceMask > 0
1110 ENDIF
1111 ENDDO
1112 ENDDO
1113
1114 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1115 #ifdef CHECK_ENERGY_CONSERV
1116 DO j = jMin, jMax
1117 DO i = iMin, iMax
1118 IF (iceMask(i,j).GT.0. _d 0) THEN
1119 qaux(1)=qIc1(i,j)
1120 qaux(2)=qIc2(i,j)
1121 CALL THSICE_CHECK_CONSERV( dBugFlag, i, j, bi, bj, 0,
1122 I iceMask(i,j), icFrac(i,j), hIce(i,j), hSnow(i,j),
1123 I qaux,
1124 I flx2oc(i,j), frw2oc(i,j), fsalt,
1125 I myTime, myIter, myThid )
1126 C endif iceMask > 0
1127 ENDIF
1128 C end i/j-loops
1129 ENDDO
1130 ENDDO
1131 #endif /* CHECK_ENERGY_CONSERV */
1132
1133 #endif /* ALLOW_THSICE */
1134
1135 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1136
1137 RETURN
1138 END

  ViewVC Help
Powered by ViewVC 1.1.22