Line data Source code
1 : !=======================================================================
2 :
3 : ! CESM meltpond parameterization
4 : !
5 : ! This meltpond parameterization was developed for use with the delta-
6 : ! Eddington radiation scheme, and only affects the radiation budget in
7 : ! the model. That is, although the pond volume is tracked, that liquid
8 : ! water is not used elsewhere in the model for mass budgets or other
9 : ! physical processes.
10 : !
11 : ! authors David A. Bailey (NCAR)
12 : ! Marika M. Holland (NCAR)
13 : ! Elizabeth C. Hunke (LANL)
14 :
15 : module icepack_meltpond_cesm
16 :
17 : use icepack_kinds
18 : use icepack_parameters, only: c0, c1, c2, p01, puny
19 : use icepack_parameters, only: rhofresh, rhoi, rhos, Timelt
20 : use icepack_warnings, only: warnstr, icepack_warnings_add
21 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
22 :
23 : implicit none
24 :
25 : private
26 : public :: compute_ponds_cesm
27 :
28 : !=======================================================================
29 :
30 : contains
31 :
32 : !=======================================================================
33 :
34 101836320 : subroutine compute_ponds_cesm(dt, hi_min, &
35 : pndaspect, &
36 : rfrac, meltt, &
37 : melts, frain, &
38 : aicen, vicen, &
39 : Tsfcn, apnd, hpnd)
40 :
41 : real (kind=dbl_kind), intent(in) :: &
42 : dt, & ! time step (s)
43 : hi_min, & ! minimum ice thickness allowed for thermo (m)
44 : pndaspect ! ratio of pond depth to pond fraction
45 :
46 : real (kind=dbl_kind), intent(in) :: &
47 : rfrac, & ! water fraction retained for melt ponds
48 : meltt, &
49 : melts, &
50 : frain, &
51 : aicen, &
52 : vicen
53 :
54 : real (kind=dbl_kind), intent(in) :: &
55 : Tsfcn
56 :
57 : real (kind=dbl_kind), intent(inout) :: &
58 : apnd, &
59 : hpnd
60 :
61 : ! local temporary variables
62 :
63 : real (kind=dbl_kind) :: &
64 960720 : volpn
65 :
66 : real (kind=dbl_kind) :: &
67 960720 : hi , & ! ice thickness (m)
68 960720 : dTs , & ! surface temperature diff for freeze-up (C)
69 960720 : Tp , & ! pond freezing temperature (C)
70 960720 : apondn, &
71 960720 : hpondn
72 :
73 : real (kind=dbl_kind), parameter :: &
74 : Td = c2 , & ! temperature difference for freeze-up (C)
75 : rexp = p01 , & ! pond contraction scaling
76 : dpthhi = 0.9_dbl_kind ! ratio of pond depth to ice thickness
77 :
78 : character(len=*),parameter :: subname='(compute_ponds_cesm)'
79 :
80 : !-----------------------------------------------------------------
81 : ! Initialize
82 : !-----------------------------------------------------------------
83 101836320 : volpn = hpnd * apnd * aicen
84 :
85 : !-----------------------------------------------------------------
86 : ! Identify grid cells where ice can melt
87 : !-----------------------------------------------------------------
88 :
89 101836320 : if (aicen > puny) then
90 :
91 22759200 : hi = vicen/aicen
92 :
93 22759200 : if (hi < hi_min) then
94 :
95 : !--------------------------------------------------------------
96 : ! Remove ponds on thin ice
97 : !--------------------------------------------------------------
98 58485 : apondn = c0
99 58485 : hpondn = c0
100 58485 : volpn = c0
101 :
102 : else
103 :
104 : !-----------------------------------------------------------
105 : ! Update pond volume
106 : !-----------------------------------------------------------
107 : volpn = volpn &
108 : + rfrac/rhofresh*(meltt*rhoi &
109 : + melts*rhos &
110 : + frain* dt)&
111 22700715 : * aicen
112 :
113 : !-----------------------------------------------------------
114 : ! Shrink pond volume under freezing conditions
115 : !-----------------------------------------------------------
116 22700715 : Tp = Timelt - Td
117 22700715 : dTs = max(Tp - Tsfcn,c0)
118 22700715 : volpn = volpn * exp(rexp*dTs/Tp)
119 22700715 : volpn = max(volpn, c0)
120 :
121 : ! fraction of ice covered by ponds
122 22700715 : apondn = min (sqrt(volpn/(pndaspect*aicen)), c1)
123 22700715 : hpondn = pndaspect * apondn
124 : ! fraction of grid cell covered by ponds
125 22700715 : apondn = apondn * aicen
126 :
127 : !-----------------------------------------------------------
128 : ! Limit pond depth
129 : !-----------------------------------------------------------
130 22700715 : hpondn = min(hpondn, dpthhi*hi)
131 :
132 : endif
133 :
134 : !-----------------------------------------------------------
135 : ! Reload tracer array
136 : !-----------------------------------------------------------
137 22759200 : apnd = apondn / aicen
138 22759200 : hpnd = hpondn
139 :
140 : endif
141 :
142 101836320 : end subroutine compute_ponds_cesm
143 :
144 : !=======================================================================
145 :
146 : end module icepack_meltpond_cesm
147 :
148 : !=======================================================================
|