Line data Source code
1 : !=======================================================================
2 :
3 : ! Flux manipulation routines for column package
4 : !
5 : ! author Elizabeth C. Hunke, LANL
6 : !
7 : ! 2014: Moved subroutines merge_fluxes, set_sfcflux from ice_flux.F90
8 :
9 : module icepack_flux
10 :
11 : use icepack_kinds
12 : use icepack_parameters, only: c1, emissivity
13 : use icepack_warnings, only: warnstr, icepack_warnings_add
14 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
15 : use icepack_tracers, only: tr_iso
16 :
17 : implicit none
18 : private
19 : public :: merge_fluxes, set_sfcflux
20 :
21 : !=======================================================================
22 :
23 : contains
24 :
25 : !=======================================================================
26 :
27 : ! Aggregate flux information from all ice thickness categories
28 : !
29 : ! author: Elizabeth C. Hunke and William H. Lipscomb, LANL
30 :
31 774714937 : subroutine merge_fluxes (aicen, &
32 : flw, &
33 : strairxn, strairyn, &
34 : Cdn_atm_ratio_n, &
35 : fsurfn, fcondtopn, &
36 : fcondbotn, &
37 : fsensn, flatn, &
38 : fswabsn, flwoutn, &
39 : evapn, &
40 : evapsn, evapin, &
41 : Trefn, Qrefn, &
42 : freshn, fsaltn, &
43 : fhocnn, fswthrun, &
44 : strairxT, strairyT, &
45 : Cdn_atm_ratio, &
46 : fsurf, fcondtop, &
47 : fcondbot, &
48 : fsens, flat, &
49 : fswabs, flwout, &
50 : evap, &
51 : evaps, evapi, &
52 : Tref, Qref, &
53 : fresh, fsalt, &
54 : fhocn, fswthru, &
55 : melttn, meltsn, meltbn, congeln, snoicen, &
56 : meltt, melts, &
57 : meltb, &
58 : congel, snoice, &
59 : Uref, Urefn, &
60 774714937 : Qref_iso, Qrefn_iso, &
61 774714937 : fiso_ocn, fiso_ocnn, &
62 774714937 : fiso_evap, fiso_evapn)
63 :
64 : ! single category fluxes
65 : real (kind=dbl_kind), intent(in) :: &
66 : aicen , & ! concentration of ice
67 : flw , & ! downward longwave flux (W/m**2)
68 : strairxn, & ! air/ice zonal strss, (N/m**2)
69 : strairyn, & ! air/ice merdnl strss, (N/m**2)
70 : Cdn_atm_ratio_n, & ! ratio of total drag over neutral drag
71 : fsurfn , & ! net heat flux to top surface (W/m**2)
72 : fcondtopn,& ! downward cond flux at top sfc (W/m**2)
73 : fcondbotn,& ! downward cond flux at bottom sfc (W/m**2)
74 : fsensn , & ! sensible heat flx (W/m**2)
75 : flatn , & ! latent heat flx (W/m**2)
76 : fswabsn , & ! shortwave absorbed heat flx (W/m**2)
77 : flwoutn , & ! upwd lw emitted heat flx (W/m**2)
78 : evapn , & ! evaporation (kg/m2/s)
79 : evapsn , & ! evaporation over snow (kg/m2/s)
80 : evapin , & ! evaporation over ice (kg/m2/s)
81 : Trefn , & ! air tmp reference level (K)
82 : Qrefn , & ! air sp hum reference level (kg/kg)
83 : freshn , & ! fresh water flux to ocean (kg/m2/s)
84 : fsaltn , & ! salt flux to ocean (kg/m2/s)
85 : fhocnn , & ! actual ocn/ice heat flx (W/m**2)
86 : fswthrun, & ! sw radiation through ice bot (W/m**2)
87 : melttn , & ! top ice melt (m)
88 : meltbn , & ! bottom ice melt (m)
89 : meltsn , & ! snow melt (m)
90 : congeln , & ! congelation ice growth (m)
91 : snoicen ! snow-ice growth (m)
92 :
93 : real (kind=dbl_kind), optional, intent(in):: &
94 : Urefn ! air speed reference level (m/s)
95 :
96 : ! cumulative fluxes
97 : real (kind=dbl_kind), intent(inout) :: &
98 : strairxT, & ! air/ice zonal strss, (N/m**2)
99 : strairyT, & ! air/ice merdnl strss, (N/m**2)
100 : Cdn_atm_ratio, & ! ratio of total drag over neutral drag
101 : fsurf , & ! net heat flux to top surface (W/m**2)
102 : fcondtop, & ! downward cond flux at top sfc (W/m**2)
103 : fcondbot, & ! downward cond flux at bottom sfc (W/m**2)
104 : fsens , & ! sensible heat flx (W/m**2)
105 : flat , & ! latent heat flx (W/m**2)
106 : fswabs , & ! shortwave absorbed heat flx (W/m**2)
107 : flwout , & ! upwd lw emitted heat flx (W/m**2)
108 : evap , & ! evaporation (kg/m2/s)
109 : evaps , & ! evaporation over snow (kg/m2/s)
110 : evapi , & ! evaporation over ice (kg/m2/s)
111 : Tref , & ! air tmp reference level (K)
112 : Qref , & ! air sp hum reference level (kg/kg)
113 : fresh , & ! fresh water flux to ocean (kg/m2/s)
114 : fsalt , & ! salt flux to ocean (kg/m2/s)
115 : fhocn , & ! actual ocn/ice heat flx (W/m**2)
116 : fswthru , & ! sw radiation through ice bot (W/m**2)
117 : meltt , & ! top ice melt (m)
118 : meltb , & ! bottom ice melt (m)
119 : melts , & ! snow melt (m)
120 : congel , & ! congelation ice growth (m)
121 : snoice ! snow-ice growth (m)
122 :
123 : real (kind=dbl_kind), optional, intent(inout):: &
124 : Uref ! air speed reference level (m/s)
125 :
126 : real (kind=dbl_kind), optional, dimension(:), intent(inout):: &
127 : Qref_iso, & ! isotope air sp hum reference level (kg/kg)
128 : fiso_ocn, & ! isotope fluxes to ocean (kg/m2/s)
129 : fiso_evap ! isotope evaporation (kg/m2/s)
130 :
131 : real (kind=dbl_kind), optional, dimension(:), intent(in):: &
132 : Qrefn_iso, & ! isotope air sp hum reference level (kg/kg)
133 : fiso_ocnn, & ! isotope fluxes to ocean (kg/m2/s)
134 : fiso_evapn ! isotope evaporation (kg/m2/s)
135 :
136 : character(len=*),parameter :: subname='(merge_fluxes)'
137 :
138 : !-----------------------------------------------------------------
139 : ! Merge fluxes
140 : ! NOTE: The albedo is aggregated only in cells where ice exists
141 : ! and (for the delta-Eddington scheme) where the sun is above
142 : ! the horizon.
143 : !-----------------------------------------------------------------
144 :
145 : ! atmo fluxes
146 :
147 774714937 : strairxT = strairxT + strairxn * aicen
148 774714937 : strairyT = strairyT + strairyn * aicen
149 : Cdn_atm_ratio = Cdn_atm_ratio + &
150 774714937 : Cdn_atm_ratio_n * aicen
151 774714937 : fsurf = fsurf + fsurfn * aicen
152 774714937 : fcondtop = fcondtop + fcondtopn * aicen
153 774714937 : fcondbot = fcondbot + fcondbotn * aicen
154 774714937 : fsens = fsens + fsensn * aicen
155 774714937 : flat = flat + flatn * aicen
156 774714937 : fswabs = fswabs + fswabsn * aicen
157 : flwout = flwout &
158 774714937 : + (flwoutn - (c1-emissivity)*flw) * aicen
159 774714937 : evap = evap + evapn * aicen
160 774714937 : evaps = evaps + evapsn * aicen
161 774714937 : evapi = evapi + evapin * aicen
162 774714937 : Tref = Tref + Trefn * aicen
163 774714937 : Qref = Qref + Qrefn * aicen
164 :
165 : ! Isotopes
166 774714937 : if (tr_iso) then
167 21958762 : if (present(Qrefn_iso) .and. present(Qref_iso)) then
168 87835048 : Qref_iso (:) = Qref_iso (:) + Qrefn_iso (:) * aicen
169 : endif
170 21958762 : if (present(fiso_ocnn) .and. present(fiso_ocn)) then
171 87835048 : fiso_ocn (:) = fiso_ocn (:) + fiso_ocnn (:) * aicen
172 : endif
173 21958762 : if (present(fiso_evapn) .and. present(fiso_evap)) then
174 87835048 : fiso_evap(:) = fiso_evap(:) + fiso_evapn(:) * aicen
175 : endif
176 : endif
177 :
178 : ! ocean fluxes
179 774714937 : if (present(Urefn) .and. present(Uref)) then
180 774714937 : Uref = Uref + Urefn * aicen
181 : endif
182 :
183 774714937 : fresh = fresh + freshn * aicen
184 774714937 : fsalt = fsalt + fsaltn * aicen
185 774714937 : fhocn = fhocn + fhocnn * aicen
186 774714937 : fswthru = fswthru + fswthrun * aicen
187 :
188 : ! ice/snow thickness
189 :
190 774714937 : meltt = meltt + melttn * aicen
191 774714937 : meltb = meltb + meltbn * aicen
192 774714937 : melts = melts + meltsn * aicen
193 774714937 : congel = congel + congeln * aicen
194 774714937 : snoice = snoice + snoicen * aicen
195 :
196 774714937 : end subroutine merge_fluxes
197 :
198 : !=======================================================================
199 :
200 : ! If model is not calculating surface temperature, set the surface
201 : ! flux values using values read in from forcing data or supplied via
202 : ! coupling (stored in ice_flux).
203 : !
204 : ! If CICE is running in NEMO environment, convert fluxes from GBM values
205 : ! to per unit ice area values. If model is not running in NEMO environment,
206 : ! the forcing is supplied as per unit ice area values.
207 : !
208 : ! authors Alison McLaren, Met Office
209 :
210 29876420 : subroutine set_sfcflux (aicen, &
211 : flatn_f, &
212 : fsensn_f, &
213 : fsurfn_f, &
214 : fcondtopn_f, &
215 : flatn, &
216 : fsensn, &
217 : fsurfn, &
218 : fcondtopn)
219 :
220 : ! ice state variables
221 : real (kind=dbl_kind), &
222 : intent(in) :: &
223 : aicen , & ! concentration of ice
224 : flatn_f , & ! latent heat flux (W/m^2)
225 : fsensn_f , & ! sensible heat flux (W/m^2)
226 : fsurfn_f , & ! net flux to top surface, not including fcondtopn
227 : fcondtopn_f ! downward cond flux at top surface (W m-2)
228 :
229 : real (kind=dbl_kind), intent(out):: &
230 : flatn , & ! latent heat flux (W/m^2)
231 : fsensn , & ! sensible heat flux (W/m^2)
232 : fsurfn , & ! net flux to top surface, not including fcondtopn
233 : fcondtopn ! downward cond flux at top surface (W m-2)
234 :
235 : ! local variables
236 :
237 : real (kind=dbl_kind) :: &
238 3773267 : raicen ! 1 or 1/aicen
239 :
240 : logical (kind=log_kind) :: &
241 : extreme_flag ! flag for extreme forcing values
242 :
243 : logical (kind=log_kind), parameter :: &
244 : extreme_test=.false. ! test and write out extreme forcing data
245 :
246 : character(len=*),parameter :: subname='(set_sfcflux)'
247 :
248 29876420 : raicen = c1
249 :
250 : #ifdef CICE_IN_NEMO
251 : !----------------------------------------------------------------------
252 : ! Convert fluxes from GBM values to per ice area values when
253 : ! running in NEMO environment. (When in standalone mode, fluxes
254 : ! are input as per ice area.)
255 : !----------------------------------------------------------------------
256 : raicen = c1 / aicen
257 : #endif
258 29876420 : fsurfn = fsurfn_f*raicen
259 29876420 : fcondtopn= fcondtopn_f*raicen
260 29876420 : flatn = flatn_f*raicen
261 29876420 : fsensn = fsensn_f*raicen
262 :
263 : !----------------------------------------------------------------
264 : ! Flag up any extreme fluxes
265 : !---------------------------------------------------------------
266 :
267 : if (extreme_test) then
268 : extreme_flag = .false.
269 :
270 : if (fcondtopn < -100.0_dbl_kind &
271 : .or. fcondtopn > 20.0_dbl_kind) then
272 : extreme_flag = .true.
273 : endif
274 :
275 : if (fsurfn < -100.0_dbl_kind &
276 : .or. fsurfn > 80.0_dbl_kind) then
277 : extreme_flag = .true.
278 : endif
279 :
280 : if (flatn < -20.0_dbl_kind &
281 : .or. flatn > 20.0_dbl_kind) then
282 : extreme_flag = .true.
283 : endif
284 :
285 : if (extreme_flag) then
286 :
287 : if (fcondtopn < -100.0_dbl_kind &
288 : .or. fcondtopn > 20.0_dbl_kind) then
289 : write(warnstr,*) subname, &
290 : 'Extreme forcing: -100 > fcondtopn > 20'
291 : call icepack_warnings_add(warnstr)
292 : write(warnstr,*) subname, &
293 : 'aicen,fcondtopn = ', &
294 : aicen,fcondtopn
295 : call icepack_warnings_add(warnstr)
296 : endif
297 :
298 : if (fsurfn < -100.0_dbl_kind &
299 : .or. fsurfn > 80.0_dbl_kind) then
300 : write(warnstr,*) subname, &
301 : 'Extreme forcing: -100 > fsurfn > 40'
302 : call icepack_warnings_add(warnstr)
303 : write(warnstr,*) subname, &
304 : 'aicen,fsurfn = ', &
305 : aicen,fsurfn
306 : call icepack_warnings_add(warnstr)
307 : endif
308 :
309 : if (flatn < -20.0_dbl_kind &
310 : .or. flatn > 20.0_dbl_kind) then
311 : write(warnstr,*) subname, &
312 : 'Extreme forcing: -20 > flatn > 20'
313 : call icepack_warnings_add(warnstr)
314 : write(warnstr,*) subname, &
315 : 'aicen,flatn = ', &
316 : aicen,flatn
317 : call icepack_warnings_add(warnstr)
318 : endif
319 :
320 : endif ! extreme_flag
321 : endif ! extreme_test
322 :
323 29876420 : end subroutine set_sfcflux
324 :
325 : !=======================================================================
326 :
327 : end module icepack_flux
328 :
329 : !=======================================================================
|