New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
limmp.F90 in branches/2016/dev_r6859_LIM3_meltponds/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2016/dev_r6859_LIM3_meltponds/NEMOGCM/NEMO/LIM_SRC_3/limmp.F90 @ 8060

Last change on this file since 8060 was 8060, checked in by vancop, 7 years ago

Simple melt ponds for interface test

File size: 87.4 KB
Line 
1MODULE limmp 
2   !!======================================================================
3   !!                     ***  MODULE  limmp   ***
4   !!   Melt ponds
5   !!======================================================================
6   !! history :       ! Original code by Daniela Flocco and Adrian Turner
7   !!            1.0  ! 2012    (O. Lecomte) Adaptation for scientific tests (NEMO3.1)
8   !!            2.0  ! 2016    (O. Lecomte, C. Rousset, M. Vancoppenolle) Implementation in NEMO3.6
9   !!----------------------------------------------------------------------
10#if defined key_lim3
11   !!----------------------------------------------------------------------
12   !!   'key_lim3' :                                 LIM3 sea-ice model
13   !!----------------------------------------------------------------------
14   !!   lim_mp_init      : some initialization and namelist read
15   !!   lim_mp           : main calling routine
16   !!   lim_mp_topo      : main melt pond routine for the "topographic" formulation (FloccoFeltham)
17   !!   lim_mp_area      : ??? compute melt pond fraction per category
18   !!   lim_mp_perm      : computes permeability (should be a FUNCTION!)
19   !!   calc_hpond       : computes melt pond depth
20   !!   permeability_phy : computes permeability
21
22   !!----------------------------------------------------------------------
23   USE phycst           ! physical constants
24   USE dom_oce          ! ocean space and time domain
25!  USE sbc_ice          ! Surface boundary condition: ice   fields
26   USE ice              ! LIM-3 variables
27   USE lbclnk           ! lateral boundary conditions - MPP exchanges
28   USE lib_mpp          ! MPP library
29   USE wrk_nemo         ! work arrays
30   USE in_out_manager   ! I/O manager
31   USE lib_fortran      ! glob_sum
32   USE timing           ! Timing
33!  USE limcons          ! conservation tests
34!  USE limctl           ! control prints
35!  USE limvar
36
37!OLI_CODE    USE ice_oce, ONLY: rdt_ice, tatm_ice
38!OLI_CODE    USE phycst
39!OLI_CODE    USE dom_ice
40!OLI_CODE    USE dom_oce
41!OLI_CODE    USE sbc_oce
42!OLI_CODE    USE sbc_ice
43!OLI_CODE    USE par_ice
44!OLI_CODE    USE par_oce
45!OLI_CODE    USE ice
46!OLI_CODE    USE thd_ice
47!OLI_CODE    USE in_out_manager
48!OLI_CODE    USE lbclnk
49!OLI_CODE    USE lib_mpp
50!OLI_CODE
51!OLI_CODE    IMPLICIT NONE
52!OLI_CODE    PRIVATE
53!OLI_CODE
54!OLI_CODE    PUBLIC   lim_mp_init
55!OLI_CODE    PUBLIC   lim_mp
56
57   IMPLICIT NONE
58   PRIVATE
59
60   PUBLIC   lim_mp_init    ! routine called by sbcice_lim.F90
61   PUBLIC   lim_mp         ! routine called by sbcice_lim.F90
62
63   !! * Substitutions
64#  include "vectopt_loop_substitute.h90"
65   !!----------------------------------------------------------------------
66   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
67   !! $Id: limdyn.F90 6994 2016-10-05 13:07:10Z clem $
68   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
69   !!----------------------------------------------------------------------
70CONTAINS
71
72   SUBROUTINE lim_mp_init 
73      !!-------------------------------------------------------------------
74      !!                  ***  ROUTINE lim_mp_init   ***
75      !!
76      !! ** Purpose : Physical constants and parameters linked to melt ponds
77      !!      over sea ice
78      !!
79      !! ** Method  :  Read the namicemp  namelist and check the melt pond 
80      !!       parameter values called at the first timestep (nit000)
81      !!
82      !! ** input   :   Namelist namicemp 
83      !!-------------------------------------------------------------------
84      INTEGER  ::   ios                 ! Local integer output status for namelist read
85      NAMELIST/namicemp/  ln_limMP, nn_limMP
86      !!-------------------------------------------------------------------
87
88      REWIND( numnam_ice_ref )              ! Namelist namicemp  in reference namelist : Melt Ponds 
89      READ  ( numnam_ice_ref, namicemp, IOSTAT = ios, ERR = 901)
90901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicemp  in reference namelist', lwp )
91
92      REWIND( numnam_ice_cfg )              ! Namelist namicemp  in configuration namelist : Melt Ponds
93      READ  ( numnam_ice_cfg, namicemp, IOSTAT = ios, ERR = 902 )
94902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicemp in configuration namelist', lwp )
95      IF(lwm) WRITE ( numoni, namicemp )
96     
97      IF(lwp) THEN                        ! control print
98         WRITE(numout,*)
99         WRITE(numout,*) 'lim_mp_init : ice parameters for melt ponds'
100         WRITE(numout,*) '~~~~~~~~~~~~'
101         WRITE(numout,*)'    Activate melt ponds                                         ln_limMP      = ', ln_limMP 
102         WRITE(numout,*)'    Type of melt pond implementation (1=all,2=rad, 3=fw, 4=no)  nn_limMP      = ', nn_limMP 
103      ENDIF
104      !
105   END SUBROUTINE lim_mp_init
106
107   SUBROUTINE lim_mp( kt )
108      !!-------------------------------------------------------------------
109      !!               ***  ROUTINE lim_mp   ***
110      !!               
111      !! ** Purpose :   change melt pond fraction
112      !!               
113      !! ** Method  :   brutal force
114      !!
115      !! ** Action  : -
116      !!              -
117      !!------------------------------------------------------------------------------------
118
119      INTEGER, INTENT(in) :: kt    ! number of iteration
120      INTEGER  ::   ji, jj, jl     ! dummy loop indices
121
122!     REAL(wp), POINTER, DIMENSION(:,:) ::   zfsurf  ! surface heat flux(obsolete, should be droped)
123      !
124      !!-------------------------------------------------------------------
125
126      IF( nn_timing == 1 )  CALL timing_start('limthd')
127
128      IF( nn_timing == 1 )  CALL timing_start('lim_mp')
129
130      CALL lim_mp_cesm ! test routine with hyper simple melt ponds
131
132      !CALL lim_mp_topo    (at_i, a_i,                                       &
133      !          &          vt_i, v_i, v_s,            t_i, s_i, a_ip_frac,  &
134      !          &          h_ip,     t_su)
135
136
137      ! we should probably not aggregate here since we do it in lim_var_agg
138      ! before output, unless we need the total volume and faction else where
139
140      ! we should also make sure a_ip and v_ip are properly updated at the end
141      ! of the routine
142
143   END SUBROUTINE lim_mp 
144
145   SUBROUTINE lim_mp_cesm
146       !!-------------------------------------------------------------------
147       !!                ***  ROUTINE lim_mp_cesm  ***
148       !!
149       !! ** Purpose    : Compute melt pond evolution
150       !!
151       !! ** Method     : Accumulation of meltwater and exponential release
152       !!
153       !! ** Tunable parameters :
154       !!               
155       !!
156       !! ** Note       : Stolen from CICE for quick test of the melt pond
157       !!                 interface
158       !!
159       !! ** References : Holland, M. M. et al (J Clim 2012)
160       !!   
161       !!-------------------------------------------------------------------
162
163       INTEGER, POINTER, DIMENSION(:,:)    :: indxi             ! compressed indices for cells with ice melting
164       INTEGER, POINTER, DIMENSION(:,:)    :: indxj             !
165
166       REAL(wp), POINTER, DIMENSION(:,:)   :: zwfx_mlw          ! available meltwater for melt ponding
167       REAL(wp), POINTER, DIMENSION(:,:,:) :: zrfrac            ! fraction of available meltwater retained for melt ponding
168
169       REAL(wp), PARAMETER                 :: zrmin  = 0.15_wp  ! minimum fraction of available meltwater retained for melt ponding
170       REAL(wp), PARAMETER                 :: zrmax  = 0.70_wp  ! maximum   ''           ''       ''        ''            ''
171       REAL(wp), PARAMETER                 :: zrexp  = 0.01_wp  ! rate constant to refreeze melt ponds
172       REAL(wp), PARAMETER                 :: zpnd_aspect = 0.8_wp ! pond aspect ratio
173
174       REAL(wp)                            :: zhi               ! dummy ice thickness
175       REAL(wp)                            :: zhs               ! dummy snow depth
176       REAL(wp)                            :: zTp               ! reference temperature
177       REAL(wp)                            :: zrexp             ! rate constant to refreeze melt ponds
178       REAL(wp)                            :: zdTs              ! dummy temperature difference
179       REAL(wp)                            :: z1_rhofw          ! inverse freshwater density
180       REAL(wp)                            :: z1_zpnd_aspect    ! inverse pond aspect ratio
181
182       !!-------------------------------------------------------------------
183
184       CALL wrk_alloc( jpi,jpj, zwfx_mlw )
185       CALL wrk_alloc( jpi*jpj, indxi, indxj)
186
187       z1_rhofw       = 1. / rhofw 
188       z1_zpnd_aspect = 1. / zpnd_aspect
189       zTp            = -2 + rt0.
190
191       !------------------------------------------------------------------
192       ! Available melt water for melt ponding and corresponding fraction
193       !------------------------------------------------------------------
194
195       zwfx_mlw(:,:) = wfx_sum(:,:) + wfx_snw(:,:)        ! available meltwater for melt ponding
196
197       zrfrac(:,:,:) = zrmin + ( zrmax - zrmin ) * a_i(:,:,:) 
198
199       DO jl = 1, jpl   
200
201          ! v_ip(:,:,jl) ! Initialize things
202          ! a_ip(:,:,jl)
203          ! volpn(:,:) = hpnd(:,:) * apnd(:,:) * aicen(:,:)
204
205          !------------------------------------------------------------------------------
206          ! Identify grid cells where ponds should be updated (can probably be improved)
207          !------------------------------------------------------------------------------
208
209          indxi(:) = 0
210          indxj(:) = 0
211          icells   = 0
212
213          DO jj = 1, jpj
214            DO ji = 1, jpi
215               IF ( a_i(ji,jj,jl) > epsi10 ) THEN
216                  icells = icells + 1
217                  indxi(icells) = i
218                  indxj(icells) = j
219               ENDIF
220            END DO                 ! ji
221         END DO                    ! jj
222
223         DO ij = 1, icells
224
225            ji = indxi(ij)
226            jj = indxj(ij)
227
228            zhi = v_i(ji,jj,jl) / a_i(ji,jj,jl)
229            zhs = v_s(ji,jj,jl) / a_i(ji,jj,jl)
230
231            IF ( zhi < rn_himin) THEN   !--- Remove ponds on thin ice if ice is too thin
232
233               wfx_pnd(ji,jj)   = wfx_pnd(ji,jj) + v_ip(ji,jj,jl) !--- Give freshwater to the ocean
234
235               a_ip(ji,jj,jl)      = 0._wp                        !--- Dump ponds
236               v_ip(ji,jj,jl)      = 0._wp
237               a_ip_frac(ji,jj,jl) = 0._wp
238               h_ip(ji,jj,jl)      = 0._wp
239
240            ELSE                        !--- Update pond characteristics
241
242               !--- Add retained melt water
243               v_ip(ji,jj,jl) = v_ip(ji,jj,jl) + rfrac(ji,jj,jl) * z1_rhofw * zwfx_mlw(ji,jj) * a_i(ji,jj,jl)
244
245               !--- Shrink pond due to refreezing
246               zdTs = MAX ( zTp - t_su(ji,jj,jl) , 0. )
247               
248               v_ip(ji,jj,jl) = v_ip(ji,jj,jl) * EXP( rexp * zdTs / zTp )
249           
250               a_ip_frac(ji,jj,jl) = MIN( 1._wp , SQRT( v_ip(ji,jj,jl) * z1_zpnd_aspect / a_i(ji,jj,jl) ) )
251
252               h_ip(ji,jj,jl)      = zpnd_aspect * a_ip_frac(ji,jj,jl)
253
254               a_ip(ji,jj,jl)      = a_ip_frac(ji,jj,jl) * a_i(ji,jj,jl)
255
256            !-----------------------------------------------------------
257            ! Limit pond depth
258            !-----------------------------------------------------------
259            ! hpondn = min(hpondn, dpthhi*hi)
260
261            !--- Give freshwater to the ocean ?
262
263            ENDIF
264
265         END DO
266
267      END DO ! jpl
268
269   END SUBROUTINE lim_mp_cesm
270
271   SUBROUTINE lim_mp_topo    (aice,      aicen,     &
272                              vice,      vicen,     &
273                              vsnon,                &
274                              ticen,     salin,     &
275                              a_ip_frac, h_ip,      &
276                                         Tsfc )
277       !!-------------------------------------------------------------------
278       !!                ***  ROUTINE lim_mp_topo  ***
279       !!
280       !! ** Purpose :   Compute melt pond evolution based on the ice
281       !!                topography as inferred from the ice thickness
282       !!                distribution. 
283       !!
284       !! ** Method  :   This code is initially based on Flocco and Feltham
285       !!                (2007) and Flocco et al. (2010). More to come...
286       !!
287       !! ** Tunable parameters :
288       !!
289       !! ** Note :
290       !!
291       !! ** References
292       !!    Flocco, D. and D. L. Feltham, 2007.  A continuum model of melt pond
293       !!    evolution on Arctic sea ice.  J. Geophys. Res. 112, C08016, doi:
294       !!    10.1029/2006JC003836.
295       !!    Flocco, D., D. L. Feltham and A. K. Turner, 2010.  Incorporation of
296       !!    a physically based melt pond scheme into the sea ice component of a
297       !!    climate model.  J. Geophys. Res. 115, C08012,
298       !!    doi: 10.1029/2009JC005568.
299       !!   
300       !!-------------------------------------------------------------------
301 
302       REAL (wp), DIMENSION (jpi,jpj), &
303          INTENT(IN) :: &
304          aice, &    ! total ice area fraction
305          vice       ! total ice volume (m)
306 
307       REAL (wp), DIMENSION (jpi,jpj,jpl), &
308          INTENT(IN) :: &
309          aicen, &   ! ice area fraction, per category
310          vsnon, &   ! snow volume, per category (m)
311          vicen      ! ice volume, per category (m)
312 
313       REAL (wp), DIMENSION (jpi,jpj,nlay_i,jpl), &
314          INTENT(IN) :: &
315          ticen, &   ! ice enthalpy, per category
316          salin
317 
318       REAL (wp), DIMENSION (jpi,jpj,jpl), &
319          INTENT(INOUT) :: &
320          a_ip_frac , &   ! pond area fraction of ice, per ice category
321          h_ip       ! pond depth, per ice category (m)
322 
323       REAL (wp), DIMENSION (jpi,jpj,jpl), &
324          INTENT(IN) :: &
325          Tsfc       ! snow/sea ice surface temperature
326 
327       ! local variables
328       REAL (wp), DIMENSION (jpi,jpj,jpl) :: &
329          zTsfcn, & ! ice/snow surface temperature (C)
330          zvolpn, & ! pond volume per unit area, per category (m)
331          zvuin     ! water-equivalent volume of ice lid on melt pond ('upper ice', m)
332 
333       REAL (wp), DIMENSION (jpi,jpj,jpl) :: &
334          zapondn,& ! pond area fraction, per category
335          zhpondn   ! pond depth, per category (m)
336 
337       REAL (wp), DIMENSION (jpi,jpj) :: &
338          zvolp       ! total volume of pond, per unit area of pond (m)
339 
340       REAL (wp) :: &
341          zhi,    & ! ice thickness (m)
342          zdHui,  & ! change in thickness of ice lid (m)
343          zomega, & ! conduction
344          zdTice, & ! temperature difference across ice lid (C)
345          zdvice, & ! change in ice volume (m)
346          zTavg,  & ! mean surface temperature across categories (C)
347          zTp,    & ! pond freezing temperature (C)
348          zdvn      ! change in melt pond volume for fresh water budget
349       INTEGER, DIMENSION (jpi*jpj) :: &
350          indxi, indxj    ! compressed indices for cells with ice melting
351 
352       INTEGER :: n,k,i,j,ij,icells,indxij ! loop indices
353 
354       INTEGER, DIMENSION (jpl) :: &
355          kcells          ! cells where ice lid combines with vice
356 
357       INTEGER, DIMENSION (jpi*jpj,jpl) :: &
358          indxii, indxjj  ! i,j indices for kcells loop
359 
360       REAL (wp), parameter :: &
361          zhicemin  = 0.1_wp , & ! minimum ice thickness with ponds (m)
362          zTd       = 0.15_wp, & ! temperature difference for freeze-up (C)
363          zr1_rlfus = 1._wp / 0.334e+6 / 917._wp , & ! (J/m^3)
364          zmin_volp = 1.e-4_wp, & ! minimum pond volume (m)
365          z0       = 0._wp,    & 
366          zTimelt   = 0._wp,    &
367          z01      = 0.01_wp,  &
368          z25      = 0.25_wp,  &
369          z5       = 0.5_wp
370
371       !---------------------------------------------------------------
372       ! Initialization
373       !---------------------------------------------------------------
374 
375       zhpondn(:,:,:) = 0._wp
376       zapondn(:,:,:) = 0._wp
377       indxii(:,:) = 0
378       indxjj(:,:) = 0
379       kcells(:)   = 0
380
381       zvolp(:,:) = wfx_sum(:,:) + wfx_snw(:,:) + vt_ip(:,:) ! Total available melt water, to be distributed as melt ponds
382       zTsfcn(:,:,:) = zTsfcn(:,:,:) - rt0                   ! Convert in Celsius
383 
384       ! The freezing temperature for meltponds is assumed slightly below 0C,
385       ! as if meltponds had a little salt in them.  The salt budget is not
386       ! altered for meltponds, but if it were then an actual pond freezing
387       ! temperature could be computed.
388 
389       ! zTp = zTimelt - zTd  ---> for lids
390 
391       !-----------------------------------------------------------------
392       ! Identify grid cells with ponds
393       !-----------------------------------------------------------------
394 
395       icells = 0
396       DO j = 1, jpj
397       DO i = 1, jpi
398          zhi = z0
399          IF (aice(i,j) > epsi10 ) zhi = vice(i,j)/aice(i,j)
400          IF ( aice(i,j) > z01 .and. zhi > zhicemin .and. &
401             zvolp(i,j) > zmin_volp*aice(i,j)) THEN
402             icells = icells + 1
403             indxi(icells) = i
404             indxj(icells) = j
405          ELSE  ! remove ponds on thin ice
406             !fpond(i,j) = fpond(i,j) - zvolp(i,j)
407             zvolpn(i,j,:) = z0
408             zvuin (i,j,:) = z0
409             zvolp (i,j) = z0
410          END IF
411       END DO                     ! i
412       END DO                     ! j
413 
414       DO ij = 1, icells
415          i = indxi(ij)
416          j = indxj(ij)
417 
418          !--------------------------------------------------------------
419          ! calculate pond area and depth
420          !--------------------------------------------------------------
421          CALL lim_mp_area(aice(i,j),vice(i,j), &
422                    aicen(i,j,:), vicen(i,j,:), vsnon(i,j,:), &
423                    ticen(i,j,:,:), salin(i,j,:,:), &
424                    zvolpn(i,j,:), zvolp(i,j), &
425                    zapondn(i,j,:),zhpondn(i,j,:), zdvn)
426         ! outputs are
427         ! - zdvn
428         ! - zvolpn
429         ! - zvolp
430         ! - zapondn
431         ! - zhpondn
432
433          wfx_pnd(i,j) = wfx_pnd(i,j) + zdvn ! update flux from ponds to ocean
434 
435          ! mean surface temperature MV - why do we need that ? --> for the lid
436
437          ! zTavg = z0
438          ! DO n = 1, jpl
439          !   zTavg = zTavg + zTsfcn(i,j,n)*aicen(i,j,n)
440          ! END DO
441          ! zTavg = zTavg / aice(i,j)
442 
443       END DO ! ij
444 
445       !---------------------------------------------------------------
446       ! Update pond volume and fraction
447       !---------------------------------------------------------------
448 
449       a_ip(:,:,:) = zapondn(:,:,:)
450       v_ip(:,:,:) = zapondn(:,:,:) * zhpondn(:,:,:)
451       a_ip_frac(:,:,:) = 0._wp
452       h_ip     (:,:,:) = 0._wp
453
454    END SUBROUTINE lim_mp_topo
455
456    SUBROUTINE lim_mp_area(aice,vice, &
457                         aicen, vicen, vsnon, ticen, &
458                         salin, zvolpn, zvolp,         &
459                         zapondn,zhpondn,dvolp)
460
461       !!-------------------------------------------------------------------
462       !!                ***  ROUTINE lim_mp_area ***
463       !!
464       !! ** Purpose : Given the total volume of meltwater, update
465       !!              pond fraction (a_ip) and depth (should be volume)
466       !!
467       !! **
468       !!             
469       !!------------------------------------------------------------------
470
471       REAL (wp), INTENT(IN) :: &
472          aice,vice
473 
474       REAL (wp), DIMENSION(jpl), INTENT(IN) :: &
475          aicen, vicen, vsnon
476 
477       REAL (wp), DIMENSION(nlay_i,jpl), INTENT(IN) :: &
478          ticen, salin
479 
480       REAL (wp), DIMENSION(jpl), INTENT(INOUT) :: &
481          zvolpn
482 
483       REAL (wp), INTENT(INOUT) :: &
484          zvolp, dvolp
485 
486       REAL (wp), DIMENSION(jpl), INTENT(OUT) :: &
487          zapondn, zhpondn
488 
489       INTEGER :: &
490          n, ns,   &
491          m_index, &
492          permflag
493 
494       REAL (wp), DIMENSION(jpl) :: &
495          hicen, &
496          hsnon, &
497          asnon, &
498          alfan, &
499          betan, &
500          cum_max_vol, &
501          reduced_aicen       
502 
503       REAL (wp), DIMENSION(0:jpl) :: &
504          cum_max_vol_tmp
505 
506       REAL (wp) :: &
507          hpond, &
508          drain, &
509          floe_weight, &
510          pressure_head, &
511          hsl_rel, &
512          deltah, &
513          perm, &
514          msno
515 
516       REAL (wp), parameter :: & 
517          viscosity = 1.79e-3_wp, &  ! kinematic water viscosity in kg/m/s
518          z0        = 0.0_wp    , &
519          c1        = 1.0_wp    , &
520          p4        = 0.4_wp    , &
521          p6        = 0.6_wp    , &
522          epsi10      = 1.0e-11_wp
523         
524      !-----------|
525      !           |
526      !           |-----------|
527      !___________|___________|______________________________________sea-level
528      !           |           |
529      !           |           |---^--------|
530      !           |           |   |        |
531      !           |           |   |        |-----------|              |-------
532      !           |           |   |alfan(n)|           |              |
533      !           |           |   |        |           |--------------|
534      !           |           |   |        |           |              |
535      !---------------------------v-------------------------------------------
536      !           |           |   ^        |           |              |
537      !           |           |   |        |           |--------------|
538      !           |           |   |betan(n)|           |              |
539      !           |           |   |        |-----------|              |-------
540      !           |           |   |        |
541      !           |           |---v------- |
542      !           |           |
543      !           |-----------|
544      !           |
545      !-----------|
546     
547       !-------------------------------------------------------------------
548       ! initialize
549       !-------------------------------------------------------------------
550 
551       DO n = 1, jpl
552 
553          zapondn(n) = z0
554          zhpondn(n) = z0
555 
556          !----------------------------------------
557          ! X) compute the effective snow fraction
558          !----------------------------------------
559          IF (aicen(n) < epsi10)  THEN
560             hicen(n) =  z0 
561             hsnon(n) = z0
562             reduced_aicen(n) = z0
563          ELSE
564             hicen(n) = vicen(n) / aicen(n)
565             hsnon(n) = vsnon(n) / aicen(n)
566             reduced_aicen(n) = c1 ! n=jpl
567             IF (n < jpl) reduced_aicen(n) = aicen(n) &
568                                  * (-0.024_wp*hicen(n) + 0.832_wp) 
569             asnon(n) = reduced_aicen(n)  ! effective snow fraction (empirical)
570             ! MV should check whether this makes sense to have the same effective snow fraction in here
571          END IF
572 
573 ! This choice for alfa and beta ignores hydrostatic equilibium of categories.
574 ! Hydrostatic equilibium of the entire ITD is accounted for below, assuming
575 ! a surface topography implied by alfa=0.6 and beta=0.4, and rigidity across all
576 ! categories.  alfa and beta partition the ITD - they are areas not thicknesses!
577 ! Multiplying by hicen, alfan and betan (below) are thus volumes per unit area.
578 ! Here, alfa = 60% of the ice area (and since hice is constant in a category,
579 ! alfan = 60% of the ice volume) in each category lies above the reference line,
580 ! and 40% below. Note: p6 is an arbitrary choice, but alfa+beta=1 is required.
581
582 ! MV: 
583 ! Note that this choice is not in the original FF07 paper and has been adopted in CICE
584 ! No reason why is explained in the doc, but I guess there is a reason. I'll try to investigate, maybe
585
586 ! Where does that choice come from
587 
588          alfan(n) = 0.6 * hicen(n)
589          betan(n) = 0.4 * hicen(n)
590       
591          cum_max_vol(n)     = z0
592          cum_max_vol_tmp(n) = z0
593     
594       END DO ! jpl
595 
596       cum_max_vol_tmp(0) = z0
597       drain = z0
598       dvolp = z0
599
600       !----------------------------------------------------------
601       ! x) Drain overflow water, update pond fraction and volume
602       !----------------------------------------------------------
603       
604       !--------------------------------------------------------------------------
605       ! the maximum amount of water that can be contained up to each ice category
606       !--------------------------------------------------------------------------
607
608       ! MV
609       ! If melt ponds are too deep to be sustainable given the ITD (OVERFLOW)
610       ! Then the excess volume cum_max_vol(jl) drains out of the system
611       ! It should be added to wfx_pnd
612       ! END MV
613     
614       DO n = 1, jpl-1 ! last category can not hold any volume
615 
616          IF (alfan(n+1) >= alfan(n) .and. alfan(n+1) > z0) THEN
617 
618             ! total volume in level including snow
619             cum_max_vol_tmp(n) = cum_max_vol_tmp(n-1) + &
620                (alfan(n+1) - alfan(n)) * sum(reduced_aicen(1:n)) 
621 
622             ! subtract snow solid volumes from lower categories in current level
623             DO ns = 1, n 
624                cum_max_vol_tmp(n) = cum_max_vol_tmp(n) &
625                   - rhosn/rhofw * &   ! free air fraction that can be filled by water
626                     asnon(ns)  * &    ! effective areal fraction of snow in that category
627                     max(min(hsnon(ns)+alfan(ns)-alfan(n), alfan(n+1)- &
628                                   alfan(n)), z0) 
629             END DO
630
631          ELSE ! assume higher categories unoccupied
632             cum_max_vol_tmp(n) = cum_max_vol_tmp(n-1)
633          END IF
634          !IF (cum_max_vol_tmp(n) < z0) THEN
635          !   call abort_ice('negative melt pond volume')
636          !END IF
637       END DO
638       cum_max_vol_tmp(jpl) = cum_max_vol_tmp(jpl-1)  ! last category holds no volume
639       cum_max_vol  (1:jpl) = cum_max_vol_tmp(1:jpl)
640     
641       !----------------------------------------------------------------
642       ! is there more meltwater than can be held in the floe?
643       !----------------------------------------------------------------
644       IF (zvolp >= cum_max_vol(jpl)) THEN
645          drain = zvolp - cum_max_vol(jpl) + epsi10
646          zvolp = zvolp - drain ! update meltwater volume available
647          dvolp = drain         ! this is the drained water
648          IF (zvolp < epsi10) THEN
649             dvolp = dvolp + zvolp
650             zvolp = z0
651          END IF
652       END IF
653     
654       ! height and area corresponding to the remaining volume
655 
656!      call calc_hpond(reduced_aicen, asnon, hsnon, rhos, alfan, &
657!           zvolp, cum_max_vol, hpond, m_index)
658     
659       DO n=1, m_index
660          zhpondn(n) = hpond - alfan(n) + alfan(1) ! here oui choulde update
661                                                   !  volume instead, no ?
662          zapondn(n) = reduced_aicen(n) 
663          ! in practise, pond fraction depends on the empirical snow fraction
664          ! so in turn on ice thickness
665       END DO
666     
667       !------------------------------------------------------------------------
668       ! Drainage through brine network (permeability)
669       !------------------------------------------------------------------------
670       !!! drainage due to ice permeability - Darcy's law
671     
672       ! sea water level
673       msno = z0
674       DO n=1,jpl
675         msno = msno + vsnon(n) * rhosn
676       END DO
677       floe_weight = (msno + rhoic*vice + rau0*zvolp) / aice
678       hsl_rel = floe_weight / rau0 &
679               - ((sum(betan(:)*aicen(:))/aice) + alfan(1))
680     
681       deltah = hpond - hsl_rel
682       pressure_head = grav * rau0 * max(deltah, z0)
683 
684       ! drain IF ice is permeable   
685       permflag = 0
686       IF (pressure_head > z0) THEN
687       DO n = 1, jpl-1
688          IF (hicen(n) /= z0) THEN
689             perm = 0. ! MV ugly dummy patch
690             CALL lim_mp_perm(ticen(:,n), salin(:,n), vicen(n), perm)
691             IF (perm > z0) permflag = 1
692
693             drain = perm*zapondn(n)*pressure_head*rdt_ice / &
694                                      (viscosity*hicen(n))
695             dvolp = dvolp + min(drain, zvolp)
696             zvolp = max(zvolp - drain, z0)
697             IF (zvolp < epsi10) THEN
698                dvolp = dvolp + zvolp
699                zvolp = z0
700             END IF
701          END IF
702       END DO
703 
704       ! adjust melt pond DIMENSIONs
705       IF (permflag > 0) THEN
706          ! recompute pond depth   
707!         CALL calc_hpond(reduced_aicen, asnon, hsnon, rhos, alfan, &
708!                         zvolp, cum_max_vol, hpond, m_index)
709          DO n=1, m_index
710             zhpondn(n) = hpond - alfan(n) + alfan(1)
711             zapondn(n) = reduced_aicen(n) 
712          END DO
713       END IF
714       END IF ! pressure_head
715 
716       !-------------------------------
717       ! X) remove water from the snow
718       !-------------------------------
719       !------------------------------------------------------------------------
720       ! total melt pond volume in category DOes not include snow volume
721       ! snow in melt ponds is not melted
722       !------------------------------------------------------------------------
723 
724       ! Calculate pond volume for lower categories
725       DO n=1,m_index-1
726          zvolpn(n) = zapondn(n) * zhpondn(n) & ! what is not in the snow
727                   - (rhosn/rhofw) * asnon(n) * min(hsnon(n), zhpondn(n))
728       END DO
729 
730       ! Calculate pond volume for highest category = remaining pond volume
731
732       ! The following is completely unclear to Martin at least
733       ! Could we redefine properly and recode in a more readable way ?
734
735       ! m_index = last category with melt pond
736
737       IF (m_index == 1) zvolpn(m_index) = zvolp ! volume of mw in 1st category is the total volume of melt water
738
739       IF (m_index > 1) THEN
740         IF (zvolp > sum(zvolpn(1:m_index-1))) THEN
741           zvolpn(m_index) = zvolp - sum(zvolpn(1:m_index-1)) !
742         ELSE
743           zvolpn(m_index) = z0
744           zhpondn(m_index) = z0
745           zapondn(m_index) = z0
746           ! If remaining pond volume is negative reduce pond volume of
747           ! lower category
748           IF (zvolp+epsi10 < sum(zvolpn(1:m_index-1))) & 
749             zvolpn(m_index-1) = zvolpn(m_index-1)-sum(zvolpn(1:m_index-1))&
750                                + zvolp
751         END IF
752       END IF
753 
754       DO n=1,m_index
755          IF (zapondn(n) > epsi10) THEN
756              zhpondn(n) = zvolpn(n) / zapondn(n)
757          ELSE
758             dvolp = dvolp + zvolpn(n)
759             zhpondn(n) = z0
760             zvolpn(n) = z0
761             zapondn(n) = z0
762          end IF
763       END DO
764       DO n = m_index+1, jpl
765          zhpondn(n) = z0
766          zapondn(n) = z0
767          zvolpn (n) = z0
768       END DO
769 
770    END SUBROUTINE lim_mp_area
771
772!OLI_CODE   
773!OLI_CODE
774!OLI_CODE    SUBROUTINE calc_hpond(aicen, asnon, hsnon, rhos, alfan, &
775!OLI_CODE                          zvolp, cum_max_vol,                &
776!OLI_CODE                          hpond, m_index)
777!OLI_CODE       !!-------------------------------------------------------------------
778!OLI_CODE       !!                ***  ROUTINE calc_hpond  ***
779!OLI_CODE       !!
780!OLI_CODE       !! ** Purpose :   Compute melt pond depth
781!OLI_CODE       !!-------------------------------------------------------------------
782!OLI_CODE     
783!OLI_CODE       REAL (wp), DIMENSION(jpl), INTENT(IN) :: &
784!OLI_CODE          aicen, &
785!OLI_CODE          asnon, &
786!OLI_CODE          hsnon, &
787!OLI_CODE          rhos,  &
788!OLI_CODE          alfan, &
789!OLI_CODE          cum_max_vol
790!OLI_CODE     
791!OLI_CODE       REAL (wp), INTENT(IN) :: &
792!OLI_CODE          zvolp
793!OLI_CODE     
794!OLI_CODE       REAL (wp), INTENT(OUT) :: &
795!OLI_CODE          hpond
796!OLI_CODE     
797!OLI_CODE       INTEGER, INTENT(OUT) :: &
798!OLI_CODE          m_index
799!OLI_CODE     
800!OLI_CODE       INTEGER :: n, ns
801!OLI_CODE     
802!OLI_CODE       REAL (wp), DIMENSION(0:jpl+1) :: &
803!OLI_CODE          hitl, &
804!OLI_CODE          aicetl
805!OLI_CODE     
806!OLI_CODE       REAL (wp) :: &
807!OLI_CODE          rem_vol, &
808!OLI_CODE          area, &
809!OLI_CODE          vol, &
810!OLI_CODE          tmp, &
811!OLI_CODE          z0   = 0.0_wp,    &   
812!OLI_CODE          epsi10 = 1.0e-11_wp
813!OLI_CODE     
814!OLI_CODE       !----------------------------------------------------------------
815!OLI_CODE       ! hpond is zero if zvolp is zero - have we fully drained?
816!OLI_CODE       !----------------------------------------------------------------
817!OLI_CODE     
818!OLI_CODE       IF (zvolp < epsi10) THEN
819!OLI_CODE        hpond = z0
820!OLI_CODE        m_index = 0
821!OLI_CODE       ELSE
822!OLI_CODE       
823!OLI_CODE        !----------------------------------------------------------------
824!OLI_CODE        ! Calculate the category where water fills up to
825!OLI_CODE        !----------------------------------------------------------------
826!OLI_CODE       
827!OLI_CODE        !----------|
828!OLI_CODE        !          |
829!OLI_CODE        !          |
830!OLI_CODE        !          |----------|                                     -- --
831!OLI_CODE        !__________|__________|_________________________________________ ^
832!OLI_CODE        !          |          |             rem_vol     ^                | Semi-filled
833!OLI_CODE        !          |          |----------|-- -- -- - ---|-- ---- -- -- --v layer
834!OLI_CODE        !          |          |          |              |
835!OLI_CODE        !          |          |          |              |hpond
836!OLI_CODE        !          |          |          |----------|   |     |-------
837!OLI_CODE        !          |          |          |          |   |     |
838!OLI_CODE        !          |          |          |          |---v-----|
839!OLI_CODE        !          |          | m_index  |          |         |
840!OLI_CODE        !-------------------------------------------------------------
841!OLI_CODE       
842!OLI_CODE        m_index = 0  ! 1:m_index categories have water in them
843!OLI_CODE        DO n = 1, jpl
844!OLI_CODE           IF (zvolp <= cum_max_vol(n)) THEN
845!OLI_CODE              m_index = n
846!OLI_CODE              IF (n == 1) THEN
847!OLI_CODE                 rem_vol = zvolp
848!OLI_CODE              ELSE
849!OLI_CODE                 rem_vol = zvolp - cum_max_vol(n-1)
850!OLI_CODE              END IF
851!OLI_CODE              exit ! to break out of the loop
852!OLI_CODE           END IF
853!OLI_CODE        END DO
854!OLI_CODE        m_index = min(jpl-1, m_index)
855!OLI_CODE       
856!OLI_CODE        !----------------------------------------------------------------
857!OLI_CODE        ! semi-filled layer may have m_index different snow in it
858!OLI_CODE        !----------------------------------------------------------------
859!OLI_CODE       
860!OLI_CODE        !-----------------------------------------------------------  ^
861!OLI_CODE        !                                                             |  alfan(m_index+1)
862!OLI_CODE        !                                                             |
863!OLI_CODE        !hitl(3)-->                             |----------|          |
864!OLI_CODE        !hitl(2)-->                |------------| * * * * *|          |
865!OLI_CODE        !hitl(1)-->     |----------|* * * * * * |* * * * * |          |
866!OLI_CODE        !hitl(0)-->-------------------------------------------------  |  ^
867!OLI_CODE        !                various snow from lower categories          |  |alfa(m_index)
868!OLI_CODE       
869!OLI_CODE        ! hitl - heights of the snow layers from thinner and current categories
870!OLI_CODE        ! aicetl - area of each snow depth in this layer
871!OLI_CODE       
872!OLI_CODE        hitl(:) = z0
873!OLI_CODE        aicetl(:) = z0
874!OLI_CODE        DO n = 1, m_index
875!OLI_CODE           hitl(n)   = max(min(hsnon(n) + alfan(n) - alfan(m_index), &
876!OLI_CODE                                  alfan(m_index+1) - alfan(m_index)), z0)
877!OLI_CODE           aicetl(n) = asnon(n)
878!OLI_CODE           
879!OLI_CODE           aicetl(0) = aicetl(0) + (aicen(n) - asnon(n))
880!OLI_CODE        END DO
881!OLI_CODE        hitl(m_index+1) = alfan(m_index+1) - alfan(m_index)
882!OLI_CODE        aicetl(m_index+1) = z0
883!OLI_CODE       
884!OLI_CODE        !----------------------------------------------------------------
885!OLI_CODE        ! reorder array according to hitl
886!OLI_CODE        ! snow heights not necessarily in height order
887!OLI_CODE        !----------------------------------------------------------------
888!OLI_CODE       
889!OLI_CODE        DO ns = 1, m_index+1
890!OLI_CODE           DO n = 0, m_index - ns + 1
891!OLI_CODE              IF (hitl(n) > hitl(n+1)) THEN ! swap order
892!OLI_CODE                 tmp = hitl(n)
893!OLI_CODE                 hitl(n) = hitl(n+1)
894!OLI_CODE                 hitl(n+1) = tmp
895!OLI_CODE                 tmp = aicetl(n)
896!OLI_CODE                 aicetl(n) = aicetl(n+1)
897!OLI_CODE                 aicetl(n+1) = tmp
898!OLI_CODE              END IF
899!OLI_CODE           END DO
900!OLI_CODE        END DO
901!OLI_CODE       
902!OLI_CODE        !----------------------------------------------------------------
903!OLI_CODE        ! divide semi-filled layer into set of sublayers each vertically homogenous
904!OLI_CODE        !----------------------------------------------------------------
905!OLI_CODE       
906!OLI_CODE        !hitl(3)----------------------------------------------------------------
907!OLI_CODE        !                                                       | * * * * * * * * 
908!OLI_CODE        !                                                       |* * * * * * * * *
909!OLI_CODE        !hitl(2)----------------------------------------------------------------
910!OLI_CODE        !                                    | * * * * * * * *  | * * * * * * * * 
911!OLI_CODE        !                                    |* * * * * * * * * |* * * * * * * * *
912!OLI_CODE        !hitl(1)----------------------------------------------------------------
913!OLI_CODE        !                 | * * * * * * * *  | * * * * * * * *  | * * * * * * * * 
914!OLI_CODE        !                 |* * * * * * * * * |* * * * * * * * * |* * * * * * * * *
915!OLI_CODE        !hitl(0)----------------------------------------------------------------
916!OLI_CODE        !    aicetl(0)         aicetl(1)           aicetl(2)          aicetl(3)           
917!OLI_CODE       
918!OLI_CODE        ! move up over layers incrementing volume
919!OLI_CODE        DO n = 1, m_index+1
920!OLI_CODE           
921!OLI_CODE           area = sum(aicetl(:)) - &                 ! total area of sub-layer
922!OLI_CODE                (rhos(n)/rau0) * sum(aicetl(n:jpl+1)) ! area of sub-layer occupied by snow
923!OLI_CODE           
924!OLI_CODE           vol = (hitl(n) - hitl(n-1)) * area      ! thickness of sub-layer times area
925!OLI_CODE           
926!OLI_CODE           IF (vol >= rem_vol) THEN  ! have reached the sub-layer with the depth within
927!OLI_CODE              hpond = rem_vol / area + hitl(n-1) + alfan(m_index) - &
928!OLI_CODE                           alfan(1)
929!OLI_CODE              exit
930!OLI_CODE           ELSE  ! still in sub-layer below the sub-layer with the depth
931!OLI_CODE              rem_vol = rem_vol - vol
932!OLI_CODE           END IF
933!OLI_CODE           
934!OLI_CODE        END DO
935!OLI_CODE       
936!OLI_CODE       END IF
937!OLI_CODE     
938!OLI_CODE    END SUBROUTINE calc_hpond
939!OLI_CODE   
940!OLI_CODE
941    SUBROUTINE lim_mp_perm(ticen, salin, vicen, perm)
942       !!-------------------------------------------------------------------
943       !!                ***  ROUTINE lim_mp_perm ***
944       !!
945       !! ** Purpose :   Determine the liquid fraction of brine in the ice
946       !!                and its permeability
947       !!-------------------------------------------------------------------
948       REAL (wp), DIMENSION(nlay_i), INTENT(IN) :: &
949          ticen,  & ! energy of melting for each ice layer (J/m2)
950          salin
951 
952       REAL (wp), INTENT(IN) :: &
953          vicen     ! ice volume
954     
955       REAL (wp), INTENT(OUT) :: &
956          perm      ! permeability
957 
958       REAL (wp) ::   &
959          Sbr        ! brine salinity
960 
961       REAL (wp), DIMENSION(nlay_i) ::   &
962          Tin, &    ! ice temperature
963          phi       ! liquid fraction
964 
965       INTEGER :: k
966     
967       REAL (wp) :: &
968          c2    = 2.0_wp
969 
970       !-----------------------------------------------------------------
971       ! Compute ice temperatures from enthalpies using quadratic formula
972       !-----------------------------------------------------------------
973 
974       DO k = 1,nlay_i
975          Tin(k) = ticen(k) 
976       END DO
977 
978       !-----------------------------------------------------------------
979       ! brine salinity and liquid fraction
980       !-----------------------------------------------------------------
981 
982       IF (maxval(Tin-rtt) <= -c2) THEN
983 
984          DO k = 1,nlay_i
985             Sbr = - 1.2_wp                 &
986                   -21.8_wp     * (Tin(k)-rtt)    &
987                   - 0.919_wp   * (Tin(k)-rtt)**2 &
988                   - 0.01878_wp * (Tin(k)-rtt)**3
989             phi(k) = salin(k)/Sbr ! liquid fraction
990          END DO ! k
991       
992       ELSE
993 
994          DO k = 1,nlay_i
995             Sbr = -17.6_wp    * (Tin(k)-rtt)    &
996                   - 0.389_wp  * (Tin(k)-rtt)**2 &
997                   - 0.00362_wp* (Tin(k)-rtt)**3
998             phi(k) = salin(k)/Sbr ! liquid fraction
999          END DO
1000 
1001       END IF
1002     
1003       !-----------------------------------------------------------------
1004       ! permeability
1005       !-----------------------------------------------------------------
1006 
1007       perm = 3.0e-08_wp * (minval(phi))**3 ! REFERENCE PLEASE (this fucking
1008                                            ! bastard of Golden)
1009     
1010    END SUBROUTINE lim_mp_perm
1011!OLI_CODE   
1012!OLI_CODE #else
1013!OLI_CODE    !!----------------------------------------------------------------------
1014!OLI_CODE    !!   Default option         Dummy Module          No LIM-3 sea-ice model
1015!OLI_CODE    !!----------------------------------------------------------------------
1016!OLI_CODE CONTAINS
1017!OLI_CODE    SUBROUTINE lim_mp_init          ! Empty routine
1018!OLI_CODE    END SUBROUTINE lim_mp_init   
1019!OLI_CODE    SUBROUTINE lim_mp               ! Empty routine
1020!OLI_CODE    END SUBROUTINE lim_mp
1021!OLI_CODE    SUBROUTINE compute_mp_topo      ! Empty routine
1022!OLI_CODE    END SUBROUTINE compute_mp_topo       
1023!OLI_CODE    SUBROUTINE pond_area            ! Empty routine
1024!OLI_CODE    END SUBROUTINE pond_area   
1025!OLI_CODE    SUBROUTINE calc_hpond           ! Empty routine
1026!OLI_CODE    END SUBROUTINE calc_hpond   
1027!OLI_CODE    SUBROUTINE permeability_phy     ! Empty routine
1028!OLI_CODE    END SUBROUTINE permeability_phy   
1029!OLI_CODE #endif
1030!OLI_CODE    !!======================================================================
1031!OLI_CODE END MODULE limmp_topo
1032
1033
1034!OLI_CODE MODULE limmp_topo
1035!OLI_CODE    !!======================================================================
1036!OLI_CODE    !!                       ***  MODULE limmp_topo ***
1037!OLI_CODE    !! LIM-3 sea-ice :  computation of melt ponds' properties
1038!OLI_CODE    !!======================================================================
1039!OLI_CODE    !! History :  Original code by Daniela Flocco and Adrian Turner
1040!OLI_CODE    !!            ! 2012-09 (O. Lecomte) Adaptation for routine inclusion in
1041!OLI_CODE    !!                      NEMO-LIM3.1
1042!OLI_CODE    !!            ! 2016-11 (O. Lecomte, C. Rousset, M. Vancoppenolle)
1043!OLI_CODE    !!                      Adaptation for merge with NEMO-LIM3.6
1044!OLI_CODE    !!----------------------------------------------------------------------
1045!OLI_CODE #if defined key_lim3
1046!OLI_CODE    !!----------------------------------------------------------------------
1047!OLI_CODE    !!   'key_lim3'                                      LIM-3 sea-ice model
1048!OLI_CODE    !!----------------------------------------------------------------------
1049!OLI_CODE    !!   lim_mp_init     : melt pond properties initialization
1050!OLI_CODE    !!   lim_mp          : melt pond routine caller
1051!OLI_CODE    !!   compute_mp_topo : Actual melt pond routine
1052!OLI_CODE    !!----------------------------------------------------------------------
1053!OLI_CODE    USE ice_oce, ONLY: rdt_ice, tatm_ice
1054!OLI_CODE    USE phycst
1055!OLI_CODE    USE dom_ice
1056!OLI_CODE    USE dom_oce
1057!OLI_CODE    USE sbc_oce
1058!OLI_CODE    USE sbc_ice
1059!OLI_CODE    USE par_ice
1060!OLI_CODE    USE par_oce
1061!OLI_CODE    USE ice
1062!OLI_CODE    USE thd_ice
1063!OLI_CODE    USE in_out_manager
1064!OLI_CODE    USE lbclnk
1065!OLI_CODE    USE lib_mpp
1066!OLI_CODE
1067!OLI_CODE    IMPLICIT NONE
1068!OLI_CODE    PRIVATE
1069!OLI_CODE
1070!OLI_CODE    PUBLIC   lim_mp_init
1071!OLI_CODE    PUBLIC   lim_mp
1072!OLI_CODE
1073!OLI_CODE CONTAINS
1074!OLI_CODE
1075!OLI_CODE    SUBROUTINE lim_mp_init
1076!OLI_CODE       !!-------------------------------------------------------------------
1077!OLI_CODE       !!                ***  ROUTINE lim_mp_init  ***
1078!OLI_CODE       !!
1079!OLI_CODE       !! ** Purpose :   Initialize melt ponds
1080!OLI_CODE       !!-------------------------------------------------------------------
1081!OLI_CODE       a_ip_frac(:,:,:)    = 0._wp
1082!OLI_CODE       a_ip(:,:,:)         = 0._wp
1083!OLI_CODE       h_ip(:,:,:)         = 0._wp
1084!OLI_CODE       v_ip(:,:,:)         = 0._wp
1085!OLI_CODE       h_il(:,:,:)         = 0._wp
1086!OLI_CODE       v_il(:,:,:)         = 0._wp
1087!OLI_CODE         
1088!OLI_CODE    END SUBROUTINE lim_mp_init
1089!OLI_CODE
1090!OLI_CODE
1091!OLI_CODE    SUBROUTINE lim_mp
1092!OLI_CODE       !!-------------------------------------------------------------------
1093!OLI_CODE       !!                ***  ROUTINE lim_mp  ***
1094!OLI_CODE       !!
1095!OLI_CODE       !! ** Purpose :   Compute surface heat flux and call main melt pond
1096!OLI_CODE       !!                routine
1097!OLI_CODE       !!-------------------------------------------------------------------
1098!OLI_CODE
1099!OLI_CODE       INTEGER  ::   ji, jj, jl   ! dummy loop indices
1100!OLI_CODE
1101!OLI_CODE       fsurf(:,:) = 0.e0
1102!OLI_CODE       DO jl = 1, jpl
1103!OLI_CODE          DO jj = 1, jpj
1104!OLI_CODE             DO ji = 1, jpi
1105!OLI_CODE                fsurf(ji,jj) = fsurf(ji,jj) + a_i(ji,jj,jl) * &
1106!OLI_CODE                      (qns_ice(ji,jj,jl) + (1.0 - izero(ji,jj,jl)) &
1107!OLI_CODE                        * qsr_ice(ji,jj,jl))
1108!OLI_CODE             END DO
1109!OLI_CODE          END DO
1110!OLI_CODE       END DO
1111!OLI_CODE
1112!OLI_CODE       CALL compute_mp_topo(at_i, a_i,                               &
1113!OLI_CODE                          vt_i, v_i, v_s, rhosn_glo, t_i, s_i, a_ip_frac,  &
1114!OLI_CODE                       h_ip, h_il, t_su, tatm_ice, diag_sur_me*rdt_ice, &
1115!OLI_CODE                          fsurf, fwoc)
1116!OLI_CODE
1117!OLI_CODE       at_ip(:,:) = 0.0
1118!OLI_CODE       vt_ip(:,:) = 0.0
1119!OLI_CODE       vt_il(:,:) = 0.0
1120!OLI_CODE       DO jl = 1, jpl
1121!OLI_CODE         DO jj = 1, jpj
1122!OLI_CODE            DO ji = 1, jpi
1123!OLI_CODE               a_ip(ji,jj,jl) = MAX(0.0_wp, a_ip_frac(ji,jj,jl) * &
1124!OLI_CODE                                                   a_i(ji,jj,jl))
1125!OLI_CODE               v_ip(ji,jj,jl) = MAX(0.0_wp, a_ip_frac(ji,jj,jl) * &
1126!OLI_CODE                                   a_i(ji,jj,jl) * h_ip(ji,jj,jl))
1127!OLI_CODE               v_il(ji,jj,jl) = MAX(0.0_wp, a_ip_frac(ji,jj,jl) * &
1128!OLI_CODE                                   a_i(ji,jj,jl) * h_il(ji,jj,jl))
1129!OLI_CODE               at_ip(ji,jj)   = at_ip(ji,jj) + a_ip(ji,jj,jl)
1130!OLI_CODE               vt_ip(ji,jj)   = vt_ip(ji,jj) + v_ip(ji,jj,jl)
1131!OLI_CODE               vt_il(ji,jj)   = vt_il(ji,jj) + v_il(ji,jj,jl)
1132!OLI_CODE            END DO
1133!OLI_CODE         END DO
1134!OLI_CODE       END DO
1135!OLI_CODE
1136!OLI_CODE    END SUBROUTINE lim_mp
1137!OLI_CODE
1138!OLI_CODE
1139!OLI_CODE    SUBROUTINE compute_mp_topo(aice,      aicen,     &
1140!OLI_CODE                               vice,      vicen,     &
1141!OLI_CODE                               vsnon,     rhos,      &
1142!OLI_CODE                               ticen,     salin,     &
1143!OLI_CODE                               a_ip_frac, h_ip,      &
1144!OLI_CODE                               h_il,      Tsfc,      &
1145!OLI_CODE                               potT,      meltt,     &
1146!OLI_CODE                               fsurf,     fwoc)
1147!OLI_CODE       !!-------------------------------------------------------------------
1148!OLI_CODE       !!                ***  ROUTINE compute_mp_topo  ***
1149!OLI_CODE       !!
1150!OLI_CODE       !! ** Purpose :   Compute melt pond evolution based on the ice
1151!OLI_CODE       !!                topography as inferred from the ice thickness
1152!OLI_CODE       !!                distribution. 
1153!OLI_CODE       !!
1154!OLI_CODE       !! ** Method  :   This code is initially based on Flocco and Feltham
1155!OLI_CODE       !!                (2007) and Flocco et al. (2010). More to come...
1156!OLI_CODE       !!
1157!OLI_CODE       !! ** Tunable parameters :
1158!OLI_CODE       !!
1159!OLI_CODE       !! ** Note :
1160!OLI_CODE       !!
1161!OLI_CODE       !! ** References
1162!OLI_CODE       !!    Flocco, D. and D. L. Feltham, 2007.  A continuum model of melt pond
1163!OLI_CODE       !!    evolution on Arctic sea ice.  J. Geophys. Res. 112, C08016, doi:
1164!OLI_CODE       !!    10.1029/2006JC003836.
1165!OLI_CODE       !!    Flocco, D., D. L. Feltham and A. K. Turner, 2010.  Incorporation of
1166!OLI_CODE       !!    a physically based melt pond scheme into the sea ice component of a
1167!OLI_CODE       !!    climate model.  J. Geophys. Res. 115, C08012,
1168!OLI_CODE       !!    doi: 10.1029/2009JC005568.
1169!OLI_CODE       !!   
1170!OLI_CODE       !!-------------------------------------------------------------------
1171!OLI_CODE
1172!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj), &
1173!OLI_CODE          INTENT(IN) :: &
1174!OLI_CODE          aice, &    ! total ice area fraction
1175!OLI_CODE          vice       ! total ice volume (m)
1176!OLI_CODE
1177!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,jpl), &
1178!OLI_CODE          INTENT(IN) :: &
1179!OLI_CODE          aicen, &   ! ice area fraction, per category
1180!OLI_CODE          vsnon, &   ! snow volume, per category (m)
1181!OLI_CODE          rhos,  &   ! equivalent snow density, per category (kg/m^3)
1182!OLI_CODE          vicen      ! ice volume, per category (m)
1183!OLI_CODE
1184!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,nlay_i,jpl), &
1185!OLI_CODE          INTENT(IN) :: &
1186!OLI_CODE          ticen, &   ! ice enthalpy, per category
1187!OLI_CODE          salin
1188!OLI_CODE
1189!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,jpl), &
1190!OLI_CODE          INTENT(INOUT) :: &
1191!OLI_CODE          a_ip_frac , &   ! pond area fraction of ice, per ice category
1192!OLI_CODE          h_ip , &   ! pond depth, per ice category (m)
1193!OLI_CODE          h_il       ! Refrozen ice lid thickness, per ice category (m)
1194!OLI_CODE
1195!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj), &
1196!OLI_CODE          INTENT(IN) :: &
1197!OLI_CODE          potT,  &   ! air potential temperature
1198!OLI_CODE          meltt, &   ! total surface meltwater flux
1199!OLI_CODE          fsurf      ! thermodynamic heat flux at ice/snow surface (W/m^2)
1200!OLI_CODE
1201!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj), &
1202!OLI_CODE          INTENT(INOUT) :: &
1203!OLI_CODE          fwoc      ! fresh water flux to the ocean (from draining and other pond volume adjustments)
1204!OLI_CODE                    ! (m)
1205!OLI_CODE
1206!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,jpl), &
1207!OLI_CODE          INTENT(IN) :: &
1208!OLI_CODE          Tsfc       ! snow/sea ice surface temperature
1209!OLI_CODE
1210!OLI_CODE       ! local variables
1211!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,jpl) :: &
1212!OLI_CODE          zTsfcn, & ! ice/snow surface temperature (C)
1213!OLI_CODE          zvolpn, & ! pond volume per unit area, per category (m)
1214!OLI_CODE          zvuin     ! water-equivalent volume of ice lid on melt pond ('upper ice', m)
1215!OLI_CODE
1216!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj,jpl) :: &
1217!OLI_CODE          zapondn,& ! pond area fraction, per category
1218!OLI_CODE          zhpondn   ! pond depth, per category (m)
1219!OLI_CODE
1220!OLI_CODE       REAL (wp), DIMENSION (jpi,jpj) :: &
1221!OLI_CODE          zvolp       ! total volume of pond, per unit area of pond (m)
1222!OLI_CODE
1223!OLI_CODE       REAL (wp) :: &
1224!OLI_CODE          zhi,    & ! ice thickness (m)
1225!OLI_CODE          zdHui,  & ! change in thickness of ice lid (m)
1226!OLI_CODE          zomega, & ! conduction
1227!OLI_CODE          zdTice, & ! temperature difference across ice lid (C)
1228!OLI_CODE          zdvice, & ! change in ice volume (m)
1229!OLI_CODE          zTavg,  & ! mean surface temperature across categories (C)
1230!OLI_CODE          zTp,    & ! pond freezing temperature (C)
1231!OLI_CODE          zdvn      ! change in melt pond volume for fresh water budget
1232!OLI_CODE       INTEGER, DIMENSION (jpi*jpj) :: &
1233!OLI_CODE          indxi, indxj    ! compressed indices for cells with ice melting
1234!OLI_CODE
1235!OLI_CODE       INTEGER :: n,k,i,j,ij,icells,indxij ! loop indices
1236!OLI_CODE
1237!OLI_CODE       INTEGER, DIMENSION (jpl) :: &
1238!OLI_CODE          kcells          ! cells where ice lid combines with vice
1239!OLI_CODE
1240!OLI_CODE       INTEGER, DIMENSION (jpi*jpj,jpl) :: &
1241!OLI_CODE          indxii, indxjj  ! i,j indices for kcells loop
1242!OLI_CODE
1243!OLI_CODE       REAL (wp), parameter :: &
1244!OLI_CODE          zhicemin  = 0.1_wp , & ! minimum ice thickness with ponds (m)
1245!OLI_CODE          zTd       = 0.15_wp, & ! temperature difference for freeze-up (C)
1246!OLI_CODE          zr1_rlfus = 1._wp / 0.334e+6 / 917._wp , & ! (J/m^3)
1247!OLI_CODE          zmin_volp = 1.e-4_wp, & ! minimum pond volume (m)
1248!OLI_CODE          z0       = 0._wp,    &
1249!OLI_CODE          zTimelt   = 0._wp,    &
1250!OLI_CODE          z01      = 0.01_wp,  &
1251!OLI_CODE          z25      = 0.25_wp,  &
1252!OLI_CODE          z5       = 0.5_wp,   &
1253!OLI_CODE          epsi10     = 1.0e-11_wp
1254!OLI_CODE       !---------------------------------------------------------------
1255!OLI_CODE       ! initialize
1256!OLI_CODE       !---------------------------------------------------------------
1257!OLI_CODE
1258!OLI_CODE       DO j = 1, jpj
1259!OLI_CODE          DO i = 1, jpi
1260!OLI_CODE             zvolp(i,j) = z0
1261!OLI_CODE          END DO
1262!OLI_CODE       END DO
1263!OLI_CODE       DO n = 1, jpl
1264!OLI_CODE          DO j = 1, jpj
1265!OLI_CODE             DO i = 1, jpi
1266!OLI_CODE                ! load tracers
1267!OLI_CODE                zvolp(i,j) = zvolp(i,j) + h_ip(i,j,n) &
1268!OLI_CODE                                      * a_ip_frac(i,j,n) * aicen(i,j,n)
1269!OLI_CODE                zTsfcn(i,j,n) = Tsfc(i,j,n) - rtt ! convert in Celsius - Oli
1270!OLI_CODE                zvuin (i,j,n) = h_il(i,j,n) &
1271!OLI_CODE                             * a_ip_frac(i,j,n) * aicen(i,j,n)
1272!OLI_CODE
1273!OLI_CODE                zhpondn(i,j,n) = z0     ! pond depth, per category
1274!OLI_CODE                zapondn(i,j,n) = z0     ! pond area,  per category
1275!OLI_CODE             END DO
1276!OLI_CODE          END DO
1277!OLI_CODE          indxii(:,n) = 0
1278!OLI_CODE          indxjj(:,n) = 0
1279!OLI_CODE          kcells  (n) = 0
1280!OLI_CODE       END DO
1281!OLI_CODE
1282!OLI_CODE       ! The freezing temperature for meltponds is assumed slightly below 0C,
1283!OLI_CODE       ! as if meltponds had a little salt in them.  The salt budget is not
1284!OLI_CODE       ! altered for meltponds, but if it were then an actual pond freezing
1285!OLI_CODE       ! temperature could be computed.
1286!OLI_CODE
1287!OLI_CODE       zTp = zTimelt - zTd
1288!OLI_CODE
1289!OLI_CODE       !-----------------------------------------------------------------
1290!OLI_CODE       ! Identify grid cells with ponds
1291!OLI_CODE       !-----------------------------------------------------------------
1292!OLI_CODE
1293!OLI_CODE       icells = 0
1294!OLI_CODE       DO j = 1, jpj
1295!OLI_CODE       DO i = 1, jpi
1296!OLI_CODE          zhi = z0
1297!OLI_CODE          IF (aice(i,j) > epsi10) zhi = vice(i,j)/aice(i,j)
1298!OLI_CODE          IF ( aice(i,j) > z01 .and. zhi > zhicemin .and. &
1299!OLI_CODE             zvolp(i,j) > zmin_volp*aice(i,j)) THEN
1300!OLI_CODE             icells = icells + 1
1301!OLI_CODE             indxi(icells) = i
1302!OLI_CODE             indxj(icells) = j
1303!OLI_CODE          ELSE  ! remove ponds on thin ice
1304!OLI_CODE             !fpond(i,j) = fpond(i,j) - zvolp(i,j)
1305!OLI_CODE             zvolpn(i,j,:) = z0
1306!OLI_CODE             zvuin (i,j,:) = z0
1307!OLI_CODE             zvolp (i,j) = z0
1308!OLI_CODE          END IF
1309!OLI_CODE       END DO                     ! i
1310!OLI_CODE       END DO                     ! j
1311!OLI_CODE
1312!OLI_CODE       DO ij = 1, icells
1313!OLI_CODE          i = indxi(ij)
1314!OLI_CODE          j = indxj(ij)
1315!OLI_CODE
1316!OLI_CODE          !--------------------------------------------------------------
1317!OLI_CODE          ! calculate pond area and depth
1318!OLI_CODE          !--------------------------------------------------------------
1319!OLI_CODE          CALL pond_area(aice(i,j),vice(i,j),rhos(i,j,:), &
1320!OLI_CODE                    aicen(i,j,:), vicen(i,j,:), vsnon(i,j,:), &
1321!OLI_CODE                    ticen(i,j,:,:), salin(i,j,:,:), &
1322!OLI_CODE                    zvolpn(i,j,:), zvolp(i,j), &
1323!OLI_CODE                    zapondn(i,j,:),zhpondn(i,j,:), zdvn)
1324!OLI_CODE
1325!OLI_CODE          fwoc(i,j) = fwoc(i,j) + zdvn ! -> Goes to fresh water budget
1326!OLI_CODE
1327!OLI_CODE          ! mean surface temperature
1328!OLI_CODE          zTavg = z0
1329!OLI_CODE          DO n = 1, jpl
1330!OLI_CODE             zTavg = zTavg + zTsfcn(i,j,n)*aicen(i,j,n)
1331!OLI_CODE          END DO
1332!OLI_CODE          zTavg = zTavg / aice(i,j)
1333!OLI_CODE
1334!OLI_CODE          DO n = 1, jpl-1
1335!OLI_CODE                       
1336!OLI_CODE             IF (zvuin(i,j,n) > epsi10) THEN
1337!OLI_CODE
1338!OLI_CODE          !----------------------------------------------------------------
1339!OLI_CODE          ! melting: floating upper ice layer melts in whole or part
1340!OLI_CODE          !----------------------------------------------------------------
1341!OLI_CODE !               IF (zTsfcn(i,j,n) > zTp) THEN
1342!OLI_CODE                IF (zTavg > zTp) THEN
1343!OLI_CODE
1344!OLI_CODE                   zdvice = min(meltt(i,j)*zapondn(i,j,n), zvuin(i,j,n))
1345!OLI_CODE                   IF (zdvice > epsi10) THEN
1346!OLI_CODE                      zvuin (i,j,n) = zvuin (i,j,n) - zdvice
1347!OLI_CODE                      zvolpn(i,j,n) = zvolpn(i,j,n) + zdvice
1348!OLI_CODE                      zvolp (i,j)   = zvolp (i,j)   + zdvice
1349!OLI_CODE                      !fwoc(i,j)   = fwoc(i,j)   + zdvice
1350!OLI_CODE                       
1351!OLI_CODE                      IF (zvuin(i,j,n) < epsi10 .and. zvolpn(i,j,n) > puny) THEN
1352!OLI_CODE                         ! ice lid melted and category is pond covered
1353!OLI_CODE                         zvolpn(i,j,n) = zvolpn(i,j,n) + zvuin(i,j,n)
1354!OLI_CODE                         !fwoc(i,j)   = fwoc(i,j)   + zvuin(i,j,n)
1355!OLI_CODE                         zvuin(i,j,n)  = z0
1356!OLI_CODE                      END IF
1357!OLI_CODE                      zhpondn(i,j,n) = zvolpn(i,j,n) / zapondn(i,j,n)
1358!OLI_CODE                   END IF
1359!OLI_CODE
1360!OLI_CODE          !----------------------------------------------------------------
1361!OLI_CODE          ! freezing: existing upper ice layer grows
1362!OLI_CODE          !----------------------------------------------------------------
1363!OLI_CODE                ELSE IF (zvolpn(i,j,n) > epsi10) THEN ! zTavg <= zTp
1364!OLI_CODE
1365!OLI_CODE                 ! dIFferential growth of base of surface floating ice layer
1366!OLI_CODE                   zdTice = max(-zTavg, z0) ! > 0
1367!OLI_CODE                   zomega = rcdic*zdTice * zr1_rlfus
1368!OLI_CODE                   zdHui = sqrt(zomega*rdt_ice + z25*(zvuin(i,j,n)/  &
1369!OLI_CODE                         aicen(i,j,n))**2)- z5 * zvuin(i,j,n)/aicen(i,j,n)
1370!OLI_CODE
1371!OLI_CODE                   zdvice = min(zdHui*zapondn(i,j,n), zvolpn(i,j,n))   
1372!OLI_CODE                   IF (zdvice > epsi10) THEN
1373!OLI_CODE                      zvuin (i,j,n) = zvuin (i,j,n) + zdvice
1374!OLI_CODE                      zvolpn(i,j,n) = zvolpn(i,j,n) - zdvice
1375!OLI_CODE                      zvolp (i,j)   = zvolp (i,j)   - zdvice
1376!OLI_CODE                      !fwoc(i,j)   = fwoc(i,j)   - zdvice
1377!OLI_CODE                      zhpondn(i,j,n) = zvolpn(i,j,n) / zapondn(i,j,n)
1378!OLI_CODE                   END IF
1379!OLI_CODE
1380!OLI_CODE                END IF ! zTavg
1381!OLI_CODE
1382!OLI_CODE          !----------------------------------------------------------------
1383!OLI_CODE          ! freezing: upper ice layer begins to form
1384!OLI_CODE          ! note: albedo does not change
1385!OLI_CODE          !----------------------------------------------------------------
1386!OLI_CODE             ELSE ! zvuin < epsi10
1387!OLI_CODE                     
1388!OLI_CODE                ! thickness of newly formed ice
1389!OLI_CODE                ! the surface temperature of a meltpond is the same as that
1390!OLI_CODE                ! of the ice underneath (0C), and the thermodynamic surface
1391!OLI_CODE                ! flux is the same
1392!OLI_CODE                zdHui = max(-fsurf(i,j)*rdt_ice*zr1_rlfus, z0)
1393!OLI_CODE                zdvice = min(zdHui*zapondn(i,j,n), zvolpn(i,j,n)) 
1394!OLI_CODE                IF (zdvice > epsi10) THEN
1395!OLI_CODE                   zvuin (i,j,n) = zdvice
1396!OLI_CODE                   zvolpn(i,j,n) = zvolpn(i,j,n) - zdvice
1397!OLI_CODE                   zvolp (i,j)   = zvolp (i,j)   - zdvice
1398!OLI_CODE                   !fwoc(i,j)   = fwoc(i,j)   - zdvice
1399!OLI_CODE                   zhpondn(i,j,n)= zvolpn(i,j,n) / zapondn(i,j,n)
1400!OLI_CODE                END IF
1401!OLI_CODE                     
1402!OLI_CODE             END IF  ! zvuin
1403!OLI_CODE
1404!OLI_CODE          END DO ! jpl
1405!OLI_CODE
1406!OLI_CODE       END DO ! ij
1407!OLI_CODE
1408!OLI_CODE       !---------------------------------------------------------------
1409!OLI_CODE       ! remove ice lid if there is no liquid pond
1410!OLI_CODE       ! zvuin may be nonzero on category jpl due to dynamics
1411!OLI_CODE       !---------------------------------------------------------------
1412!OLI_CODE
1413!OLI_CODE       DO j = 1, jpj
1414!OLI_CODE       DO i = 1, jpi
1415!OLI_CODE          DO n = 1, jpl
1416!OLI_CODE             IF (aicen(i,j,n) > epsi10 .and. zvolpn(i,j,n) < puny &
1417!OLI_CODE                                     .and. zvuin (i,j,n) > epsi10) THEN
1418!OLI_CODE                kcells(n) = kcells(n) + 1
1419!OLI_CODE                indxij    = kcells(n)
1420!OLI_CODE                indxii(indxij,n) = i
1421!OLI_CODE                indxjj(indxij,n) = j
1422!OLI_CODE             END IF
1423!OLI_CODE          END DO
1424!OLI_CODE       END DO                     ! i
1425!OLI_CODE       END DO                     ! j
1426!OLI_CODE
1427!OLI_CODE       DO n = 1, jpl
1428!OLI_CODE
1429!OLI_CODE          IF (kcells(n) > 0) THEN
1430!OLI_CODE          DO ij = 1, kcells(n)
1431!OLI_CODE             i = indxii(ij,n)
1432!OLI_CODE             j = indxjj(ij,n)
1433!OLI_CODE             fwoc(i,j) = fwoc(i,j) + rhoic/rauw * zvuin(i,j,n) ! Completely refrozen lid goes into ocean (to be changed)
1434!OLI_CODE             zvuin(i,j,n) = z0
1435!OLI_CODE          END DO    ! ij
1436!OLI_CODE          END IF
1437!OLI_CODE
1438!OLI_CODE          ! reload tracers
1439!OLI_CODE          DO j = 1, jpj
1440!OLI_CODE             DO i = 1, jpi
1441!OLI_CODE                IF (zapondn(i,j,n) > epsi10) THEN
1442!OLI_CODE                   h_il(i,j,n) = zvuin(i,j,n) / zapondn(i,j,n)
1443!OLI_CODE                ELSE
1444!OLI_CODE                   zvuin(i,j,n) = z0
1445!OLI_CODE                   h_il(i,j,n) = z0
1446!OLI_CODE                END IF
1447!OLI_CODE                IF (aicen(i,j,n) > epsi10) THEN
1448!OLI_CODE                   a_ip_frac(i,j,n) = zapondn(i,j,n) / aicen(i,j,n) * &
1449!OLI_CODE                         (1.0_wp - MAX(z0, SIGN(1.0_wp, -zvolpn(i,j,n))))
1450!OLI_CODE                   h_ip(i,j,n) = zhpondn(i,j,n)
1451!OLI_CODE                ELSE
1452!OLI_CODE                   a_ip_frac(i,j,n) = z0
1453!OLI_CODE                   h_ip(i,j,n) = z0
1454!OLI_CODE                   h_il(i,j,n) = z0
1455!OLI_CODE                END IF
1456!OLI_CODE             END DO ! i
1457!OLI_CODE          END DO    ! j
1458!OLI_CODE
1459!OLI_CODE       END DO       ! n
1460!OLI_CODE
1461!OLI_CODE    END SUBROUTINE compute_mp_topo
1462!OLI_CODE
1463!OLI_CODE
1464!OLI_CODE    SUBROUTINE pond_area(aice,vice,rhos,             &
1465!OLI_CODE                         aicen, vicen, vsnon, ticen, &
1466!OLI_CODE                         salin, zvolpn, zvolp,         &
1467!OLI_CODE                         zapondn,zhpondn,dvolp)
1468!OLI_CODE       !!-------------------------------------------------------------------
1469!OLI_CODE       !!                ***  ROUTINE pond_area  ***
1470!OLI_CODE       !!
1471!OLI_CODE       !! ** Purpose :   Compute melt pond area, depth and melting rates
1472!OLI_CODE       !!------------------------------------------------------------------
1473!OLI_CODE       REAL (wp), INTENT(IN) :: &
1474!OLI_CODE          aice,vice
1475!OLI_CODE
1476!OLI_CODE       REAL (wp), DIMENSION(jpl), INTENT(IN) :: &
1477!OLI_CODE          aicen, vicen, vsnon, rhos
1478!OLI_CODE
1479!OLI_CODE       REAL (wp), DIMENSION(nlay_i,jpl), INTENT(IN) :: &
1480!OLI_CODE          ticen, salin
1481!OLI_CODE
1482!OLI_CODE       REAL (wp), DIMENSION(jpl), INTENT(INOUT) :: &
1483!OLI_CODE          zvolpn
1484!OLI_CODE
1485!OLI_CODE       REAL (wp), INTENT(INOUT) :: &
1486!OLI_CODE          zvolp, dvolp
1487!OLI_CODE
1488!OLI_CODE       REAL (wp), DIMENSION(jpl), INTENT(OUT) :: &
1489!OLI_CODE          zapondn, zhpondn
1490!OLI_CODE
1491!OLI_CODE       INTEGER :: &
1492!OLI_CODE          n, ns,   &
1493!OLI_CODE          m_index, &
1494!OLI_CODE          permflag
1495!OLI_CODE
1496!OLI_CODE       REAL (wp), DIMENSION(jpl) :: &
1497!OLI_CODE          hicen, &
1498!OLI_CODE          hsnon, &
1499!OLI_CODE          asnon, &
1500!OLI_CODE          alfan, &
1501!OLI_CODE          betan, &
1502!OLI_CODE          cum_max_vol, &
1503!OLI_CODE          reduced_aicen       
1504!OLI_CODE
1505!OLI_CODE       REAL (wp), DIMENSION(0:jpl) :: &
1506!OLI_CODE          cum_max_vol_tmp
1507!OLI_CODE
1508!OLI_CODE       REAL (wp) :: &
1509!OLI_CODE          hpond, &
1510!OLI_CODE          drain, &
1511!OLI_CODE          floe_weight, &
1512!OLI_CODE          pressure_head, &
1513!OLI_CODE          hsl_rel, &
1514!OLI_CODE          deltah, &
1515!OLI_CODE          perm, &
1516!OLI_CODE          apond, &
1517!OLI_CODE          msno
1518!OLI_CODE
1519!OLI_CODE       REAL (wp), parameter :: &
1520!OLI_CODE          viscosity = 1.79e-3_wp, &  ! kinematic water viscosity in kg/m/s
1521!OLI_CODE          z0        = 0.0_wp    , &
1522!OLI_CODE          c1        = 1.0_wp    , &
1523!OLI_CODE          p4        = 0.4_wp    , &
1524!OLI_CODE          p6        = 0.6_wp    , &
1525!OLI_CODE          epsi10      = 1.0e-11_wp
1526!OLI_CODE         
1527!OLI_CODE      !-----------|
1528!OLI_CODE      !           |
1529!OLI_CODE      !           |-----------|
1530!OLI_CODE      !___________|___________|______________________________________sea-level
1531!OLI_CODE      !           |           |
1532!OLI_CODE      !           |           |---^--------|
1533!OLI_CODE      !           |           |   |        |
1534!OLI_CODE      !           |           |   |        |-----------|              |-------
1535!OLI_CODE      !           |           |   |alfan(n)|           |              |
1536!OLI_CODE      !           |           |   |        |           |--------------|
1537!OLI_CODE      !           |           |   |        |           |              |
1538!OLI_CODE      !---------------------------v-------------------------------------------
1539!OLI_CODE      !           |           |   ^        |           |              |
1540!OLI_CODE      !           |           |   |        |           |--------------|
1541!OLI_CODE      !           |           |   |betan(n)|           |              |
1542!OLI_CODE      !           |           |   |        |-----------|              |-------
1543!OLI_CODE      !           |           |   |        |
1544!OLI_CODE      !           |           |---v------- |
1545!OLI_CODE      !           |           |
1546!OLI_CODE      !           |-----------|
1547!OLI_CODE      !           |
1548!OLI_CODE      !-----------|
1549!OLI_CODE     
1550!OLI_CODE       !-------------------------------------------------------------------
1551!OLI_CODE       ! initialize
1552!OLI_CODE       !-------------------------------------------------------------------
1553!OLI_CODE
1554!OLI_CODE       DO n = 1, jpl
1555!OLI_CODE
1556!OLI_CODE          zapondn(n) = z0
1557!OLI_CODE          zhpondn(n) = z0
1558!OLI_CODE
1559!OLI_CODE          IF (aicen(n) < epsi10)  THEN
1560!OLI_CODE             hicen(n) =  z0
1561!OLI_CODE             hsnon(n) = z0
1562!OLI_CODE             reduced_aicen(n) = z0
1563!OLI_CODE          ELSE
1564!OLI_CODE             hicen(n) = vicen(n) / aicen(n)
1565!OLI_CODE             hsnon(n) = vsnon(n) / aicen(n)
1566!OLI_CODE             reduced_aicen(n) = c1 ! n=jpl
1567!OLI_CODE             IF (n < jpl) reduced_aicen(n) = aicen(n) &
1568!OLI_CODE                                  * (-0.024_wp*hicen(n) + 0.832_wp)
1569!OLI_CODE             asnon(n) = reduced_aicen(n)
1570!OLI_CODE          END IF
1571!OLI_CODE
1572!OLI_CODE ! This choice for alfa and beta ignores hydrostatic equilibium of categories.
1573!OLI_CODE ! Hydrostatic equilibium of the entire ITD is accounted for below, assuming
1574!OLI_CODE ! a surface topography implied by alfa=0.6 and beta=0.4, and rigidity across all
1575!OLI_CODE ! categories.  alfa and beta partition the ITD - they are areas not thicknesses!
1576!OLI_CODE ! Multiplying by hicen, alfan and betan (below) are thus volumes per unit area.
1577!OLI_CODE ! Here, alfa = 60% of the ice area (and since hice is constant in a category,
1578!OLI_CODE ! alfan = 60% of the ice volume) in each category lies above the reference line,
1579!OLI_CODE ! and 40% below. Note: p6 is an arbitrary choice, but alfa+beta=1 is required.
1580!OLI_CODE
1581!OLI_CODE          alfan(n) = p6 * hicen(n)
1582!OLI_CODE          betan(n) = p4 * hicen(n)
1583!OLI_CODE       
1584!OLI_CODE          cum_max_vol(n)     = z0
1585!OLI_CODE          cum_max_vol_tmp(n) = z0
1586!OLI_CODE     
1587!OLI_CODE       END DO ! jpl
1588!OLI_CODE
1589!OLI_CODE       cum_max_vol_tmp(0) = z0
1590!OLI_CODE       drain = z0
1591!OLI_CODE       dvolp = z0
1592!OLI_CODE     
1593!OLI_CODE       !--------------------------------------------------------------------------
1594!OLI_CODE       ! the maximum amount of water that can be contained up to each ice category
1595!OLI_CODE       !--------------------------------------------------------------------------
1596!OLI_CODE     
1597!OLI_CODE       DO n = 1, jpl-1 ! last category can not hold any volume
1598!OLI_CODE
1599!OLI_CODE          IF (alfan(n+1) >= alfan(n) .and. alfan(n+1) > z0) THEN
1600!OLI_CODE
1601!OLI_CODE             ! total volume in level including snow
1602!OLI_CODE             cum_max_vol_tmp(n) = cum_max_vol_tmp(n-1) + &
1603!OLI_CODE                (alfan(n+1) - alfan(n)) * sum(reduced_aicen(1:n))
1604!OLI_CODE
1605!OLI_CODE
1606!OLI_CODE             ! subtract snow solid volumes from lower categories in current level
1607!OLI_CODE             DO ns = 1, n
1608!OLI_CODE                cum_max_vol_tmp(n) = cum_max_vol_tmp(n) &
1609!OLI_CODE                   - rhos(ns)/rauw  * &    ! fraction of snow that is occupied by solid ??rauw
1610!OLI_CODE                     asnon(ns)  * &    ! area of snow from that category
1611!OLI_CODE                     max(min(hsnon(ns)+alfan(ns)-alfan(n), alfan(n+1)- &
1612!OLI_CODE                                   alfan(n)), z0) 
1613!OLI_CODE                                       ! thickness of snow from ns layer in n layer
1614!OLI_CODE             END DO
1615!OLI_CODE
1616!OLI_CODE          ELSE ! assume higher categories unoccupied
1617!OLI_CODE             cum_max_vol_tmp(n) = cum_max_vol_tmp(n-1)
1618!OLI_CODE          END IF
1619!OLI_CODE          !IF (cum_max_vol_tmp(n) < z0) THEN
1620!OLI_CODE          !   call abort_ice('negative melt pond volume')
1621!OLI_CODE          !END IF
1622!OLI_CODE       END DO
1623!OLI_CODE       cum_max_vol_tmp(jpl) = cum_max_vol_tmp(jpl-1)  ! last category holds no volume
1624!OLI_CODE       cum_max_vol  (1:jpl) = cum_max_vol_tmp(1:jpl)
1625!OLI_CODE     
1626!OLI_CODE       !----------------------------------------------------------------
1627!OLI_CODE       ! is there more meltwater than can be held in the floe?
1628!OLI_CODE       !----------------------------------------------------------------
1629!OLI_CODE       IF (zvolp >= cum_max_vol(jpl)) THEN
1630!OLI_CODE          drain = zvolp - cum_max_vol(jpl) + epsi10
1631!OLI_CODE          zvolp = zvolp - drain
1632!OLI_CODE          dvolp = drain
1633!OLI_CODE          IF (zvolp < epsi10) THEN
1634!OLI_CODE             dvolp = dvolp + zvolp
1635!OLI_CODE             zvolp = z0
1636!OLI_CODE          END IF
1637!OLI_CODE       END IF
1638!OLI_CODE     
1639!OLI_CODE       ! height and area corresponding to the remaining volume
1640!OLI_CODE
1641!OLI_CODE       call calc_hpond(reduced_aicen, asnon, hsnon, rhos, alfan, &
1642!OLI_CODE            zvolp, cum_max_vol, hpond, m_index)
1643!OLI_CODE     
1644!OLI_CODE       DO n=1, m_index
1645!OLI_CODE          zhpondn(n) = hpond - alfan(n) + alfan(1)
1646!OLI_CODE          zapondn(n) = reduced_aicen(n)
1647!OLI_CODE       END DO
1648!OLI_CODE       apond = sum(zapondn(1:m_index))
1649!OLI_CODE     
1650!OLI_CODE       !------------------------------------------------------------------------
1651!OLI_CODE       ! drainage due to ice permeability - Darcy's law
1652!OLI_CODE       !------------------------------------------------------------------------
1653!OLI_CODE     
1654!OLI_CODE       ! sea water level
1655!OLI_CODE       msno = z0
1656!OLI_CODE       DO n=1,jpl
1657!OLI_CODE         msno = msno + vsnon(n) * rhos(n)
1658!OLI_CODE       END DO
1659!OLI_CODE       floe_weight = (msno + rhoic*vice + rau0*zvolp) / aice
1660!OLI_CODE       hsl_rel = floe_weight / rau0 &
1661!OLI_CODE               - ((sum(betan(:)*aicen(:))/aice) + alfan(1))
1662!OLI_CODE     
1663!OLI_CODE       deltah = hpond - hsl_rel
1664!OLI_CODE       pressure_head = grav * rau0 * max(deltah, z0)
1665!OLI_CODE
1666!OLI_CODE       ! drain IF ice is permeable   
1667!OLI_CODE       permflag = 0
1668!OLI_CODE       IF (pressure_head > z0) THEN
1669!OLI_CODE       DO n = 1, jpl-1
1670!OLI_CODE          IF (hicen(n) /= z0) THEN
1671!OLI_CODE             CALL permeability_phi(ticen(:,n), salin(:,n), vicen(n), perm)
1672!OLI_CODE             IF (perm > z0) permflag = 1
1673!OLI_CODE             drain = perm*zapondn(n)*pressure_head*rdt_ice / &
1674!OLI_CODE                                      (viscosity*hicen(n))
1675!OLI_CODE             dvolp = dvolp + min(drain, zvolp)
1676!OLI_CODE             zvolp = max(zvolp - drain, z0)
1677!OLI_CODE             IF (zvolp < epsi10) THEN
1678!OLI_CODE                dvolp = dvolp + zvolp
1679!OLI_CODE                zvolp = z0
1680!OLI_CODE             END IF
1681!OLI_CODE          END IF
1682!OLI_CODE       END DO
1683!OLI_CODE 
1684!OLI_CODE       ! adjust melt pond DIMENSIONs
1685!OLI_CODE       IF (permflag > 0) THEN
1686!OLI_CODE          ! recompute pond depth   
1687!OLI_CODE          CALL calc_hpond(reduced_aicen, asnon, hsnon, rhos, alfan, &
1688!OLI_CODE                          zvolp, cum_max_vol, hpond, m_index)
1689!OLI_CODE          DO n=1, m_index
1690!OLI_CODE             zhpondn(n) = hpond - alfan(n) + alfan(1)
1691!OLI_CODE             zapondn(n) = reduced_aicen(n)
1692!OLI_CODE          END DO
1693!OLI_CODE          apond = sum(zapondn(1:m_index))
1694!OLI_CODE       END IF
1695!OLI_CODE       END IF ! pressure_head
1696!OLI_CODE
1697!OLI_CODE       !------------------------------------------------------------------------
1698!OLI_CODE       ! total melt pond volume in category DOes not include snow volume
1699!OLI_CODE       ! snow in melt ponds is not melted
1700!OLI_CODE       !------------------------------------------------------------------------
1701!OLI_CODE
1702!OLI_CODE       ! Calculate pond volume for lower categories
1703!OLI_CODE       DO n=1,m_index-1
1704!OLI_CODE          zvolpn(n) = zapondn(n) * zhpondn(n) &
1705!OLI_CODE                   - (rhos(n)/rauw) * asnon(n) * min(hsnon(n), zhpondn(n))!
1706!OLI_CODE       END DO
1707!OLI_CODE
1708!OLI_CODE       ! Calculate pond volume for highest category = remaining pond volume
1709!OLI_CODE       IF (m_index == 1) zvolpn(m_index) = zvolp
1710!OLI_CODE       IF (m_index > 1) THEN
1711!OLI_CODE         IF (zvolp > sum(zvolpn(1:m_index-1))) THEN
1712!OLI_CODE           zvolpn(m_index) = zvolp - sum(zvolpn(1:m_index-1))
1713!OLI_CODE         ELSE
1714!OLI_CODE           zvolpn(m_index) = z0
1715!OLI_CODE           zhpondn(m_index) = z0
1716!OLI_CODE           zapondn(m_index) = z0
1717!OLI_CODE           ! If remaining pond volume is negative reduce pond volume of
1718!OLI_CODE           ! lower category
1719!OLI_CODE           IF (zvolp+epsi10 < sum(zvolpn(1:m_index-1))) &
1720!OLI_CODE             zvolpn(m_index-1) = zvolpn(m_index-1)-sum(zvolpn(1:m_index-1))&
1721!OLI_CODE                                + zvolp
1722!OLI_CODE         END IF
1723!OLI_CODE       END IF
1724!OLI_CODE
1725!OLI_CODE       DO n=1,m_index
1726!OLI_CODE          IF (zapondn(n) > epsi10) THEN
1727!OLI_CODE              zhpondn(n) = zvolpn(n) / zapondn(n)
1728!OLI_CODE          ELSE
1729!OLI_CODE             dvolp = dvolp + zvolpn(n)
1730!OLI_CODE             zhpondn(n) = z0
1731!OLI_CODE             zvolpn(n) = z0
1732!OLI_CODE             zapondn(n) = z0
1733!OLI_CODE          end IF
1734!OLI_CODE       END DO
1735!OLI_CODE       DO n = m_index+1, jpl
1736!OLI_CODE          zhpondn(n) = z0
1737!OLI_CODE          zapondn(n) = z0
1738!OLI_CODE          zvolpn (n) = z0
1739!OLI_CODE       END DO
1740!OLI_CODE
1741!OLI_CODE    END SUBROUTINE pond_area
1742!OLI_CODE   
1743!OLI_CODE
1744!OLI_CODE    SUBROUTINE calc_hpond(aicen, asnon, hsnon, rhos, alfan, &
1745!OLI_CODE                          zvolp, cum_max_vol,                &
1746!OLI_CODE                          hpond, m_index)
1747!OLI_CODE       !!-------------------------------------------------------------------
1748!OLI_CODE       !!                ***  ROUTINE calc_hpond  ***
1749!OLI_CODE       !!
1750!OLI_CODE       !! ** Purpose :   Compute melt pond depth
1751!OLI_CODE       !!-------------------------------------------------------------------
1752!OLI_CODE     
1753!OLI_CODE       REAL (wp), DIMENSION(jpl), INTENT(IN) :: &
1754!OLI_CODE          aicen, &
1755!OLI_CODE          asnon, &
1756!OLI_CODE          hsnon, &
1757!OLI_CODE          rhos,  &
1758!OLI_CODE          alfan, &
1759!OLI_CODE          cum_max_vol
1760!OLI_CODE     
1761!OLI_CODE       REAL (wp), INTENT(IN) :: &
1762!OLI_CODE          zvolp
1763!OLI_CODE     
1764!OLI_CODE       REAL (wp), INTENT(OUT) :: &
1765!OLI_CODE          hpond
1766!OLI_CODE     
1767!OLI_CODE       INTEGER, INTENT(OUT) :: &
1768!OLI_CODE          m_index
1769!OLI_CODE     
1770!OLI_CODE       INTEGER :: n, ns
1771!OLI_CODE     
1772!OLI_CODE       REAL (wp), DIMENSION(0:jpl+1) :: &
1773!OLI_CODE          hitl, &
1774!OLI_CODE          aicetl
1775!OLI_CODE     
1776!OLI_CODE       REAL (wp) :: &
1777!OLI_CODE          rem_vol, &
1778!OLI_CODE          area, &
1779!OLI_CODE          vol, &
1780!OLI_CODE          tmp, &
1781!OLI_CODE          z0   = 0.0_wp,    &   
1782!OLI_CODE          epsi10 = 1.0e-11_wp
1783!OLI_CODE     
1784!OLI_CODE       !----------------------------------------------------------------
1785!OLI_CODE       ! hpond is zero if zvolp is zero - have we fully drained?
1786!OLI_CODE       !----------------------------------------------------------------
1787!OLI_CODE     
1788!OLI_CODE       IF (zvolp < epsi10) THEN
1789!OLI_CODE        hpond = z0
1790!OLI_CODE        m_index = 0
1791!OLI_CODE       ELSE
1792!OLI_CODE       
1793!OLI_CODE        !----------------------------------------------------------------
1794!OLI_CODE        ! Calculate the category where water fills up to
1795!OLI_CODE        !----------------------------------------------------------------
1796!OLI_CODE       
1797!OLI_CODE        !----------|
1798!OLI_CODE        !          |
1799!OLI_CODE        !          |
1800!OLI_CODE        !          |----------|                                     -- --
1801!OLI_CODE        !__________|__________|_________________________________________ ^
1802!OLI_CODE        !          |          |             rem_vol     ^                | Semi-filled
1803!OLI_CODE        !          |          |----------|-- -- -- - ---|-- ---- -- -- --v layer
1804!OLI_CODE        !          |          |          |              |
1805!OLI_CODE        !          |          |          |              |hpond
1806!OLI_CODE        !          |          |          |----------|   |     |-------
1807!OLI_CODE        !          |          |          |          |   |     |
1808!OLI_CODE        !          |          |          |          |---v-----|
1809!OLI_CODE        !          |          | m_index  |          |         |
1810!OLI_CODE        !-------------------------------------------------------------
1811!OLI_CODE       
1812!OLI_CODE        m_index = 0  ! 1:m_index categories have water in them
1813!OLI_CODE        DO n = 1, jpl
1814!OLI_CODE           IF (zvolp <= cum_max_vol(n)) THEN
1815!OLI_CODE              m_index = n
1816!OLI_CODE              IF (n == 1) THEN
1817!OLI_CODE                 rem_vol = zvolp
1818!OLI_CODE              ELSE
1819!OLI_CODE                 rem_vol = zvolp - cum_max_vol(n-1)
1820!OLI_CODE              END IF
1821!OLI_CODE              exit ! to break out of the loop
1822!OLI_CODE           END IF
1823!OLI_CODE        END DO
1824!OLI_CODE        m_index = min(jpl-1, m_index)
1825!OLI_CODE       
1826!OLI_CODE        !----------------------------------------------------------------
1827!OLI_CODE        ! semi-filled layer may have m_index different snow in it
1828!OLI_CODE        !----------------------------------------------------------------
1829!OLI_CODE       
1830!OLI_CODE        !-----------------------------------------------------------  ^
1831!OLI_CODE        !                                                             |  alfan(m_index+1)
1832!OLI_CODE        !                                                             |
1833!OLI_CODE        !hitl(3)-->                             |----------|          |
1834!OLI_CODE        !hitl(2)-->                |------------| * * * * *|          |
1835!OLI_CODE        !hitl(1)-->     |----------|* * * * * * |* * * * * |          |
1836!OLI_CODE        !hitl(0)-->-------------------------------------------------  |  ^
1837!OLI_CODE        !                various snow from lower categories          |  |alfa(m_index)
1838!OLI_CODE       
1839!OLI_CODE        ! hitl - heights of the snow layers from thinner and current categories
1840!OLI_CODE        ! aicetl - area of each snow depth in this layer
1841!OLI_CODE       
1842!OLI_CODE        hitl(:) = z0
1843!OLI_CODE        aicetl(:) = z0
1844!OLI_CODE        DO n = 1, m_index
1845!OLI_CODE           hitl(n)   = max(min(hsnon(n) + alfan(n) - alfan(m_index), &
1846!OLI_CODE                                  alfan(m_index+1) - alfan(m_index)), z0)
1847!OLI_CODE           aicetl(n) = asnon(n)
1848!OLI_CODE           
1849!OLI_CODE           aicetl(0) = aicetl(0) + (aicen(n) - asnon(n))
1850!OLI_CODE        END DO
1851!OLI_CODE        hitl(m_index+1) = alfan(m_index+1) - alfan(m_index)
1852!OLI_CODE        aicetl(m_index+1) = z0
1853!OLI_CODE       
1854!OLI_CODE        !----------------------------------------------------------------
1855!OLI_CODE        ! reorder array according to hitl
1856!OLI_CODE        ! snow heights not necessarily in height order
1857!OLI_CODE        !----------------------------------------------------------------
1858!OLI_CODE       
1859!OLI_CODE        DO ns = 1, m_index+1
1860!OLI_CODE           DO n = 0, m_index - ns + 1
1861!OLI_CODE              IF (hitl(n) > hitl(n+1)) THEN ! swap order
1862!OLI_CODE                 tmp = hitl(n)
1863!OLI_CODE                 hitl(n) = hitl(n+1)
1864!OLI_CODE                 hitl(n+1) = tmp
1865!OLI_CODE                 tmp = aicetl(n)
1866!OLI_CODE                 aicetl(n) = aicetl(n+1)
1867!OLI_CODE                 aicetl(n+1) = tmp
1868!OLI_CODE              END IF
1869!OLI_CODE           END DO
1870!OLI_CODE        END DO
1871!OLI_CODE       
1872!OLI_CODE        !----------------------------------------------------------------
1873!OLI_CODE        ! divide semi-filled layer into set of sublayers each vertically homogenous
1874!OLI_CODE        !----------------------------------------------------------------
1875!OLI_CODE       
1876!OLI_CODE        !hitl(3)----------------------------------------------------------------
1877!OLI_CODE        !                                                       | * * * * * * * * 
1878!OLI_CODE        !                                                       |* * * * * * * * *
1879!OLI_CODE        !hitl(2)----------------------------------------------------------------
1880!OLI_CODE        !                                    | * * * * * * * *  | * * * * * * * * 
1881!OLI_CODE        !                                    |* * * * * * * * * |* * * * * * * * *
1882!OLI_CODE        !hitl(1)----------------------------------------------------------------
1883!OLI_CODE        !                 | * * * * * * * *  | * * * * * * * *  | * * * * * * * * 
1884!OLI_CODE        !                 |* * * * * * * * * |* * * * * * * * * |* * * * * * * * *
1885!OLI_CODE        !hitl(0)----------------------------------------------------------------
1886!OLI_CODE        !    aicetl(0)         aicetl(1)           aicetl(2)          aicetl(3)           
1887!OLI_CODE       
1888!OLI_CODE        ! move up over layers incrementing volume
1889!OLI_CODE        DO n = 1, m_index+1
1890!OLI_CODE           
1891!OLI_CODE           area = sum(aicetl(:)) - &                 ! total area of sub-layer
1892!OLI_CODE                (rhos(n)/rau0) * sum(aicetl(n:jpl+1)) ! area of sub-layer occupied by snow
1893!OLI_CODE           
1894!OLI_CODE           vol = (hitl(n) - hitl(n-1)) * area      ! thickness of sub-layer times area
1895!OLI_CODE           
1896!OLI_CODE           IF (vol >= rem_vol) THEN  ! have reached the sub-layer with the depth within
1897!OLI_CODE              hpond = rem_vol / area + hitl(n-1) + alfan(m_index) - &
1898!OLI_CODE                           alfan(1)
1899!OLI_CODE              exit
1900!OLI_CODE           ELSE  ! still in sub-layer below the sub-layer with the depth
1901!OLI_CODE              rem_vol = rem_vol - vol
1902!OLI_CODE           END IF
1903!OLI_CODE           
1904!OLI_CODE        END DO
1905!OLI_CODE       
1906!OLI_CODE       END IF
1907!OLI_CODE     
1908!OLI_CODE    END SUBROUTINE calc_hpond
1909!OLI_CODE   
1910!OLI_CODE
1911!OLI_CODE    SUBROUTINE permeability_phi(ticen, salin, vicen, perm)
1912!OLI_CODE       !!-------------------------------------------------------------------
1913!OLI_CODE       !!                ***  ROUTINE permeability_phi  ***
1914!OLI_CODE       !!
1915!OLI_CODE       !! ** Purpose :   Determine the liquid fraction of brine in the ice
1916!OLI_CODE       !!                and its permeability
1917!OLI_CODE       !!-------------------------------------------------------------------
1918!OLI_CODE       REAL (wp), DIMENSION(nlay_i), INTENT(IN) :: &
1919!OLI_CODE          ticen,  & ! energy of melting for each ice layer (J/m2)
1920!OLI_CODE          salin
1921!OLI_CODE
1922!OLI_CODE       REAL (wp), INTENT(IN) :: &
1923!OLI_CODE          vicen     ! ice volume
1924!OLI_CODE     
1925!OLI_CODE       REAL (wp), INTENT(OUT) :: &
1926!OLI_CODE          perm      ! permeability
1927!OLI_CODE
1928!OLI_CODE       REAL (wp) ::   &
1929!OLI_CODE          Sbr        ! brine salinity
1930!OLI_CODE
1931!OLI_CODE       REAL (wp), DIMENSION(nlay_i) ::   &
1932!OLI_CODE          Tin, &    ! ice temperature
1933!OLI_CODE          phi       ! liquid fraction
1934!OLI_CODE
1935!OLI_CODE       INTEGER :: k
1936!OLI_CODE     
1937!OLI_CODE       REAL (wp) :: &
1938!OLI_CODE          c2    = 2.0_wp
1939!OLI_CODE 
1940!OLI_CODE       !-----------------------------------------------------------------
1941!OLI_CODE       ! Compute ice temperatures from enthalpies using quadratic formula
1942!OLI_CODE       !-----------------------------------------------------------------
1943!OLI_CODE
1944!OLI_CODE       DO k = 1,nlay_i
1945!OLI_CODE          Tin(k) = ticen(k)
1946!OLI_CODE       END DO
1947!OLI_CODE
1948!OLI_CODE       !-----------------------------------------------------------------
1949!OLI_CODE       ! brine salinity and liquid fraction
1950!OLI_CODE       !-----------------------------------------------------------------
1951!OLI_CODE
1952!OLI_CODE       IF (maxval(Tin-rtt) <= -c2) THEN
1953!OLI_CODE
1954!OLI_CODE          DO k = 1,nlay_i
1955!OLI_CODE             Sbr = - 1.2_wp                 &
1956!OLI_CODE                   -21.8_wp     * (Tin(k)-rtt)    &
1957!OLI_CODE                   - 0.919_wp   * (Tin(k)-rtt)**2 &
1958!OLI_CODE                   - 0.01878_wp * (Tin(k)-rtt)**3
1959!OLI_CODE             phi(k) = salin(k)/Sbr ! liquid fraction
1960!OLI_CODE          END DO ! k
1961!OLI_CODE       
1962!OLI_CODE       ELSE
1963!OLI_CODE
1964!OLI_CODE          DO k = 1,nlay_i
1965!OLI_CODE             Sbr = -17.6_wp    * (Tin(k)-rtt)    &
1966!OLI_CODE                   - 0.389_wp  * (Tin(k)-rtt)**2 &
1967!OLI_CODE                   - 0.00362_wp* (Tin(k)-rtt)**3
1968!OLI_CODE             phi(k) = salin(k)/Sbr ! liquid fraction
1969!OLI_CODE          END DO
1970!OLI_CODE
1971!OLI_CODE       END IF
1972!OLI_CODE     
1973!OLI_CODE       !-----------------------------------------------------------------
1974!OLI_CODE       ! permeability
1975!OLI_CODE       !-----------------------------------------------------------------
1976!OLI_CODE
1977!OLI_CODE       perm = 3.0e-08_wp * (minval(phi))**3
1978!OLI_CODE     
1979!OLI_CODE    END SUBROUTINE permeability_phi
1980!OLI_CODE   
1981!OLI_CODE #else
1982!OLI_CODE    !!----------------------------------------------------------------------
1983!OLI_CODE    !!   Default option         Dummy Module          No LIM-3 sea-ice model
1984!OLI_CODE    !!----------------------------------------------------------------------
1985!OLI_CODE CONTAINS
1986!OLI_CODE    SUBROUTINE lim_mp_init          ! Empty routine
1987!OLI_CODE    END SUBROUTINE lim_mp_init   
1988!OLI_CODE    SUBROUTINE lim_mp               ! Empty routine
1989!OLI_CODE    END SUBROUTINE lim_mp
1990!OLI_CODE    SUBROUTINE compute_mp_topo      ! Empty routine
1991!OLI_CODE    END SUBROUTINE compute_mp_topo       
1992!OLI_CODE    SUBROUTINE pond_area            ! Empty routine
1993!OLI_CODE    END SUBROUTINE pond_area   
1994!OLI_CODE    SUBROUTINE calc_hpond           ! Empty routine
1995!OLI_CODE    END SUBROUTINE calc_hpond   
1996!OLI_CODE    SUBROUTINE permeability_phy     ! Empty routine
1997!OLI_CODE    END SUBROUTINE permeability_phy   
1998!OLI_CODE #endif
1999!OLI_CODE    !!======================================================================
2000!OLI_CODE END MODULE limmp_topo
2001!OLI_CODE
2002#else
2003   !!----------------------------------------------------------------------
2004   !!   Default option          Empty module           NO LIM sea-ice model
2005   !!----------------------------------------------------------------------
2006CONTAINS
2007   SUBROUTINE lim_mp_init     ! Empty routine
2008   END SUBROUTINE lim_mp_init
2009   SUBROUTINE lim_mp          ! Empty routine
2010   END SUBROUTINE lim_mp     
2011   SUBROUTINE lim_mp_topo     ! Empty routine
2012   END SUBROUTINE lim_mp_topo
2013   SUBROUTINE lim_mp_cesm     ! Empty routine
2014   END SUBROUTINE lim_mp_cesm
2015   SUBROUTINE lim_mp_area     ! Empty routine
2016   END SUBROUTINE lim_mp_area
2017   SUBROUTINE lim_mp_perm     ! Empty routine
2018   END SUBROUTINE lim_mp_perm
2019#endif 
2020
2021   !!======================================================================
2022END MODULE limmp 
Note: See TracBrowser for help on using the repository browser.