Skip to content
Snippets Groups Projects
icethd_pnd.F90 65.8 KiB
Newer Older
Guillaume Samson's avatar
Guillaume Samson committed
1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410

        !------------------------------------------------------------------------
        ! Drainage through brine network (permeability)
        !------------------------------------------------------------------------
        !!! drainage due to ice permeability - Darcy's law

        ! sea water level
        msno = 0._wp
        DO jl = 1 , jpl
          msno = msno + v_s(ji,jj,jl) * rhos
        END DO
        floe_weight = ( msno + rhoi*vt_i(ji,jj) + rho0*zvolp(ji,jj) ) / at_i(ji,jj)
        hsl_rel = floe_weight / rho0 &
                - ( ( sum(betan(:)*a_i(ji,jj,:)) / at_i(ji,jj) ) + alfan(1) )

        deltah = hpond - hsl_rel
        pressure_head = grav * rho0 * max(deltah, 0._wp)

        ! drain if ice is permeable
        permflag = 0

        IF (pressure_head > 0._wp) THEN
           DO jl = 1, jpl-1
              IF ( hicen(jl) /= 0._wp ) THEN

              !IF (hicen(jl) > 0._wp) THEN           !js: from CICE 5.1.2

                 perm = 0._wp ! MV ugly dummy patch
                 CALL ice_thd_pnd_perm(t_i(ji,jj,:,jl),  sz_i(ji,jj,:,jl), perm) ! bof
                 IF (perm > 0._wp) permflag = 1

                 drain = perm*a_ip(ji,jj,jl)*pressure_head*rDt_ice / &
                                          (viscosity*hicen(jl))
                 zdvolp(ji,jj) = zdvolp(ji,jj) + min(drain, zvolp(ji,jj))
                 zvolp(ji,jj) = max(zvolp(ji,jj) - drain, 0._wp)

                 diag_dvpn_drn(ji,jj) = - drain ! diag (could be better coded)

                 IF (zvolp(ji,jj) < epsi10) THEN
                    zdvolp(ji,jj) = zdvolp(ji,jj) + zvolp(ji,jj)
                    zvolp(ji,jj) = 0._wp
                 END IF
             END IF
          END DO

           ! adjust melt pond dimensions
           IF (permflag > 0) THEN
              ! recompute pond depth
             CALL ice_thd_pnd_depth(reduced_aicen, asnon, hsnon, alfan, zvolp(ji,jj), cum_max_vol, hpond, m_index)
              DO jl = 1, m_index
                 h_ip(ji,jj,jl) = hpond - alfan(jl) + alfan(1)
                 a_ip(ji,jj,jl) = reduced_aicen(jl)
              END DO
              !zapond = sum(a_ip(1:m_index))       !js: from CICE 5.1.2; not in Icepack1.1.0-6-gac6195d
           END IF
        END IF ! pressure_head

        !-------------------------------
        ! remove water from the snow
        !-------------------------------
        !------------------------------------------------------------------------
        ! total melt pond volume in category does not include snow volume
        ! snow in melt ponds is not melted
        !------------------------------------------------------------------------

        ! MV here, it seems that we remove some meltwater from the ponds, but I can't really tell
        ! how much, so I did not diagnose it
        ! so if there is a problem here, nobody is going to see it...


        ! Calculate pond volume for lower categories
        DO jl = 1,m_index-1
           v_ip(ji,jj,jl) = a_ip(ji,jj,jl) * h_ip(ji,jj,jl) & ! what is not in the snow
                          - (rhos/rhow) * asnon(jl) * min(hsnon(jl), h_ip(ji,jj,jl))
        END DO

        ! Calculate pond volume for highest category = remaining pond volume

        ! The following is completely unclear to Martin at least
        ! Could we redefine properly and recode in a more readable way ?

        ! m_index = last category with melt pond

        IF (m_index == 1) v_ip(ji,jj,m_index) = zvolp(ji,jj) ! volume of mw in 1st category is the total volume of melt water

        IF (m_index > 1) THEN
          IF (zvolp(ji,jj) > sum( v_ip(ji,jj,1:m_index-1))) THEN
             v_ip(ji,jj,m_index) = zvolp(ji,jj) - sum(v_ip(ji,jj,1:m_index-1))
          ELSE
             v_ip(ji,jj,m_index) = 0._wp
             h_ip(ji,jj,m_index) = 0._wp
             a_ip(ji,jj,m_index) = 0._wp
             ! If remaining pond volume is negative reduce pond volume of
             ! lower category
             IF ( zvolp(ji,jj) + epsi10 < SUM(v_ip(ji,jj,1:m_index-1))) &
              v_ip(ji,jj,m_index-1) = v_ip(ji,jj,m_index-1) - sum(v_ip(ji,jj,1:m_index-1)) + zvolp(ji,jj)
          END IF
        END IF

        DO jl = 1,m_index
           IF (a_ip(ji,jj,jl) > epsi10) THEN
               h_ip(ji,jj,jl) = v_ip(ji,jj,jl) / a_ip(ji,jj,jl)
           ELSE
              zdvolp(ji,jj) = zdvolp(ji,jj) + v_ip(ji,jj,jl)
              h_ip(ji,jj,jl) = 0._wp
              v_ip(ji,jj,jl)  = 0._wp
              a_ip(ji,jj,jl) = 0._wp
           END IF
        END DO
        DO jl = m_index+1, jpl
           h_ip(ji,jj,jl) = 0._wp
           a_ip(ji,jj,jl) = 0._wp
           v_ip(ji,jj,jl) = 0._wp
        END DO

           ENDIF

     END_2D

    END SUBROUTINE ice_thd_pnd_area


    SUBROUTINE ice_thd_pnd_depth(aicen, asnon, hsnon, alfan, zvolp, cum_max_vol, hpond, m_index)
       !!-------------------------------------------------------------------
       !!                ***  ROUTINE ice_thd_pnd_depth  ***
       !!
       !! ** Purpose :   Compute melt pond depth
       !!-------------------------------------------------------------------

       REAL (wp), DIMENSION(jpl), INTENT(IN) :: &
          aicen, &
          asnon, &
          hsnon, &
          alfan, &
          cum_max_vol

       REAL (wp), INTENT(IN) :: &
          zvolp

       REAL (wp), INTENT(OUT) :: &
          hpond

       INTEGER, INTENT(OUT) :: &
          m_index

       INTEGER :: n, ns

       REAL (wp), DIMENSION(0:jpl+1) :: &
          hitl, &
          aicetl

       REAL (wp) :: &
          rem_vol, &
          area, &
          vol, &
          tmp, &
          z0   = 0.0_wp

       !----------------------------------------------------------------
       ! hpond is zero if zvolp is zero - have we fully drained?
       !----------------------------------------------------------------

       IF (zvolp < epsi10) THEN
        hpond = z0
        m_index = 0
       ELSE

        !----------------------------------------------------------------
        ! Calculate the category where water fills up to
        !----------------------------------------------------------------

        !----------|
        !          |
        !          |
        !          |----------|                                     -- --
        !__________|__________|_________________________________________ ^
        !          |          |             rem_vol     ^                | Semi-filled
        !          |          |----------|-- -- -- - ---|-- ---- -- -- --v layer
        !          |          |          |              |
        !          |          |          |              |hpond
        !          |          |          |----------|   |     |-------
        !          |          |          |          |   |     |
        !          |          |          |          |---v-----|
        !          |          | m_index  |          |         |
        !-------------------------------------------------------------

        m_index = 0  ! 1:m_index categories have water in them
        DO n = 1, jpl
           IF (zvolp <= cum_max_vol(n)) THEN
              m_index = n
              IF (n == 1) THEN
                 rem_vol = zvolp
              ELSE
                 rem_vol = zvolp - cum_max_vol(n-1)
              END IF
              exit ! to break out of the loop
           END IF
        END DO
        m_index = min(jpl-1, m_index)

        !----------------------------------------------------------------
        ! semi-filled layer may have m_index different snow in it
        !----------------------------------------------------------------

        !-----------------------------------------------------------  ^
        !                                                             |  alfan(m_index+1)
        !                                                             |
        !hitl(3)-->                             |----------|          |
        !hitl(2)-->                |------------| * * * * *|          |
        !hitl(1)-->     |----------|* * * * * * |* * * * * |          |
        !hitl(0)-->-------------------------------------------------  |  ^
        !                various snow from lower categories          |  |alfa(m_index)

        ! hitl - heights of the snow layers from thinner and current categories
        ! aicetl - area of each snow depth in this layer

        hitl(:) = z0
        aicetl(:) = z0
        DO n = 1, m_index
           hitl(n)   = max(min(hsnon(n) + alfan(n) - alfan(m_index), &
                                  alfan(m_index+1) - alfan(m_index)), z0)
           aicetl(n) = asnon(n)

           aicetl(0) = aicetl(0) + (aicen(n) - asnon(n))
        END DO

        hitl(m_index+1) = alfan(m_index+1) - alfan(m_index)
        aicetl(m_index+1) = z0

        !----------------------------------------------------------------
        ! reorder array according to hitl
        ! snow heights not necessarily in height order
        !----------------------------------------------------------------

        DO ns = 1, m_index+1
           DO n = 0, m_index - ns + 1
              IF (hitl(n) > hitl(n+1)) THEN ! swap order
                 tmp = hitl(n)
                 hitl(n) = hitl(n+1)
                 hitl(n+1) = tmp
                 tmp = aicetl(n)
                 aicetl(n) = aicetl(n+1)
                 aicetl(n+1) = tmp
              END IF
           END DO
        END DO

        !----------------------------------------------------------------
        ! divide semi-filled layer into set of sublayers each vertically homogenous
        !----------------------------------------------------------------

        !hitl(3)----------------------------------------------------------------
        !                                                       | * * * * * * * *
        !                                                       |* * * * * * * * *
        !hitl(2)----------------------------------------------------------------
        !                                    | * * * * * * * *  | * * * * * * * *
        !                                    |* * * * * * * * * |* * * * * * * * *
        !hitl(1)----------------------------------------------------------------
        !                 | * * * * * * * *  | * * * * * * * *  | * * * * * * * *
        !                 |* * * * * * * * * |* * * * * * * * * |* * * * * * * * *
        !hitl(0)----------------------------------------------------------------
        !    aicetl(0)         aicetl(1)           aicetl(2)          aicetl(3)

        ! move up over layers incrementing volume
        DO n = 1, m_index+1

           area = sum(aicetl(:)) - &                 ! total area of sub-layer
                (rhos/rho0) * sum(aicetl(n:jpl+1)) ! area of sub-layer occupied by snow

           vol = (hitl(n) - hitl(n-1)) * area      ! thickness of sub-layer times area

           IF (vol >= rem_vol) THEN  ! have reached the sub-layer with the depth within
              hpond = rem_vol / area + hitl(n-1) + alfan(m_index) - alfan(1)

              exit
           ELSE  ! still in sub-layer below the sub-layer with the depth
              rem_vol = rem_vol - vol
           END IF

        END DO

       END IF

    END SUBROUTINE ice_thd_pnd_depth


    SUBROUTINE ice_thd_pnd_perm(ticen, salin, perm)
       !!-------------------------------------------------------------------
       !!                ***  ROUTINE ice_thd_pnd_perm ***
       !!
       !! ** Purpose :   Determine the liquid fraction of brine in the ice
       !!                and its permeability
       !!-------------------------------------------------------------------

       REAL (wp), DIMENSION(nlay_i), INTENT(IN) :: &
          ticen,  &     ! internal ice temperature (K)
          salin         ! salinity (ppt)     !js: ppt according to cice

       REAL (wp), INTENT(OUT) :: &
          perm      ! permeability

       REAL (wp) ::   &
          Sbr       ! brine salinity

       REAL (wp), DIMENSION(nlay_i) ::   &
          Tin, &    ! ice temperature
          phi       ! liquid fraction

       INTEGER :: k

       !-----------------------------------------------------------------
       ! Compute ice temperatures from enthalpies using quadratic formula
       !-----------------------------------------------------------------

       DO k = 1,nlay_i
          Tin(k) = ticen(k) - rt0   !js: from K to degC
       END DO

       !-----------------------------------------------------------------
       ! brine salinity and liquid fraction
       !-----------------------------------------------------------------

       DO k = 1, nlay_i

          Sbr    = - Tin(k) / rTmlt ! Consistent expression with SI3 (linear liquidus)
          ! Best expression to date is that one (Vancoppenolle et al JGR 2019)
          ! Sbr  = - 18.7 * Tin(k) - 0.519 * Tin(k)**2 - 0.00535 * Tin(k) **3
          phi(k) = salin(k) / Sbr

       END DO

       !-----------------------------------------------------------------
       ! permeability
       !-----------------------------------------------------------------

       perm = 3.0e-08_wp * (minval(phi))**3 ! Golden et al. (2007)

   END SUBROUTINE ice_thd_pnd_perm

   SUBROUTINE ice_thd_pnd_init
      !!-------------------------------------------------------------------
      !!                  ***  ROUTINE ice_thd_pnd_init   ***
      !!
      !! ** Purpose : Physical constants and parameters linked to melt ponds
      !!              over sea ice
      !!
      !! ** Method  :  Read the namthd_pnd  namelist and check the melt pond
      !!               parameter values called at the first timestep (nit000)
      !!
      !! ** input   :   Namelist namthd_pnd
      !!-------------------------------------------------------------------
      INTEGER  ::   ios, ioptio   ! Local integer
      !!
      NAMELIST/namthd_pnd/  ln_pnd, ln_pnd_LEV , rn_apnd_min, rn_apnd_max, rn_pnd_flush, &
         &                          ln_pnd_CST , rn_apnd, rn_hpnd,         &
         &                          ln_pnd_TOPO,                           &
         &                          ln_pnd_lids, ln_pnd_alb
      !!-------------------------------------------------------------------
      !
      READ  ( numnam_ice_ref, namthd_pnd, IOSTAT = ios, ERR = 901)
901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namthd_pnd  in reference namelist' )
      READ  ( numnam_ice_cfg, namthd_pnd, IOSTAT = ios, ERR = 902 )
902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namthd_pnd in configuration namelist' )
      IF(lwm) WRITE ( numoni, namthd_pnd )
      !
      IF(lwp) THEN                        ! control print
         WRITE(numout,*)
         WRITE(numout,*) 'ice_thd_pnd_init: ice parameters for melt ponds'
         WRITE(numout,*) '~~~~~~~~~~~~~~~~'
         WRITE(numout,*) '   Namelist namicethd_pnd:'
         WRITE(numout,*) '      Melt ponds activated or not                                 ln_pnd       = ', ln_pnd
         WRITE(numout,*) '         Topographic melt pond scheme                             ln_pnd_TOPO  = ', ln_pnd_TOPO
         WRITE(numout,*) '         Level ice melt pond scheme                               ln_pnd_LEV   = ', ln_pnd_LEV
         WRITE(numout,*) '            Minimum ice fraction that contributes to melt ponds   rn_apnd_min  = ', rn_apnd_min
         WRITE(numout,*) '            Maximum ice fraction that contributes to melt ponds   rn_apnd_max  = ', rn_apnd_max
         WRITE(numout,*) '            Pond flushing efficiency                              rn_pnd_flush = ', rn_pnd_flush
         WRITE(numout,*) '         Constant ice melt pond scheme                            ln_pnd_CST   = ', ln_pnd_CST
         WRITE(numout,*) '            Prescribed pond fraction                              rn_apnd      = ', rn_apnd
         WRITE(numout,*) '            Prescribed pond depth                                 rn_hpnd      = ', rn_hpnd
         WRITE(numout,*) '         Frozen lids on top of melt ponds                         ln_pnd_lids  = ', ln_pnd_lids
         WRITE(numout,*) '         Melt ponds affect albedo or not                          ln_pnd_alb   = ', ln_pnd_alb
      ENDIF
      !
      !                             !== set the choice of ice pond scheme ==!
      ioptio = 0
      IF( .NOT.ln_pnd ) THEN   ;   ioptio = ioptio + 1   ;   nice_pnd = np_pndNO     ;   ENDIF
      IF( ln_pnd_CST  ) THEN   ;   ioptio = ioptio + 1   ;   nice_pnd = np_pndCST    ;   ENDIF
      IF( ln_pnd_LEV  ) THEN   ;   ioptio = ioptio + 1   ;   nice_pnd = np_pndLEV    ;   ENDIF
      IF( ln_pnd_TOPO ) THEN   ;   ioptio = ioptio + 1   ;   nice_pnd = np_pndTOPO   ;   ENDIF
      IF( ioptio /= 1 )   &
         & CALL ctl_stop( 'ice_thd_pnd_init: choose either none (ln_pnd=F) or only one pond scheme (ln_pnd_LEV, ln_pnd_CST or ln_pnd_TOPO)' )
      !
      SELECT CASE( nice_pnd )
      CASE( np_pndNO )
         IF( ln_pnd_alb  ) THEN ; ln_pnd_alb  = .FALSE. ; CALL ctl_warn( 'ln_pnd_alb=false when no ponds' )  ; ENDIF
         IF( ln_pnd_lids ) THEN ; ln_pnd_lids = .FALSE. ; CALL ctl_warn( 'ln_pnd_lids=false when no ponds' ) ; ENDIF
      CASE( np_pndCST )
         IF( ln_pnd_lids ) THEN ; ln_pnd_lids = .FALSE. ; CALL ctl_warn( 'ln_pnd_lids=false when constant ponds' ) ; ENDIF
      END SELECT
      !
   END SUBROUTINE ice_thd_pnd_init

#else
   !!----------------------------------------------------------------------
   !!   Default option          Empty module          NO SI3 sea-ice model
   !!----------------------------------------------------------------------
#endif

   !!======================================================================
END MODULE icethd_pnd