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, snwgrain
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 8777938 : 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 : fswthrun_vdr, fswthrun_vdf,&
45 : fswthrun_idr, fswthrun_idf,&
46 : strairxT, strairyT, &
47 : Cdn_atm_ratio, &
48 : fsurf, fcondtop, &
49 : fcondbot, &
50 : fsens, flat, &
51 : fswabs, flwout, &
52 : evap, &
53 : evaps, evapi, &
54 : Tref, Qref, &
55 : fresh, fsalt, &
56 : fhocn, fswthru, &
57 : fswthru_vdr, fswthru_vdf,&
58 : fswthru_idr, fswthru_idf,&
59 : melttn, meltsn, meltbn, congeln, snoicen, &
60 : meltt, melts, &
61 : meltb, dsnow, dsnown,&
62 : congel, snoice, &
63 : meltsliq, meltsliqn, &
64 : Uref, Urefn, &
65 8777938 : Qref_iso, Qrefn_iso, &
66 8777938 : fiso_ocn, fiso_ocnn, &
67 8777938 : fiso_evap, fiso_evapn)
68 :
69 : ! single category fluxes
70 : real (kind=dbl_kind), intent(in) :: &
71 : aicen ! concentration of ice
72 :
73 : real (kind=dbl_kind), optional, intent(in) :: &
74 : flw , & ! downward longwave flux (W/m**2)
75 : strairxn, & ! air/ice zonal strss, (N/m**2)
76 : strairyn, & ! air/ice merdnl strss, (N/m**2)
77 : Cdn_atm_ratio_n, & ! ratio of total drag over neutral drag
78 : fsurfn , & ! net heat flux to top surface (W/m**2)
79 : fcondtopn,& ! downward cond flux at top sfc (W/m**2)
80 : fcondbotn,& ! downward cond flux at bottom sfc (W/m**2)
81 : fsensn , & ! sensible heat flx (W/m**2)
82 : flatn , & ! latent heat flx (W/m**2)
83 : fswabsn , & ! shortwave absorbed heat flx (W/m**2)
84 : flwoutn , & ! upwd lw emitted heat flx (W/m**2)
85 : evapn , & ! evaporation (kg/m2/s)
86 : evapsn , & ! evaporation over snow (kg/m2/s)
87 : evapin , & ! evaporation over ice (kg/m2/s)
88 : Trefn , & ! air tmp reference level (K)
89 : Qrefn , & ! air sp hum reference level (kg/kg)
90 : freshn , & ! fresh water flux to ocean (kg/m2/s)
91 : fsaltn , & ! salt flux to ocean (kg/m2/s)
92 : fhocnn , & ! actual ocn/ice heat flx (W/m**2)
93 : fswthrun, & ! sw radiation through ice bot (W/m**2)
94 : melttn , & ! top ice melt (m)
95 : meltbn , & ! bottom ice melt (m)
96 : meltsn , & ! snow melt (m)
97 : meltsliqn,& ! mass of snow melt (kg/m^2)
98 : dsnown , & ! change in snow depth (m)
99 : congeln , & ! congelation ice growth (m)
100 : snoicen , & ! snow-ice growth (m)
101 : fswthrun_vdr, & ! vis dir sw radiation through ice bot (W/m**2)
102 : fswthrun_vdf, & ! vis dif sw radiation through ice bot (W/m**2)
103 : fswthrun_idr, & ! nir dir sw radiation through ice bot (W/m**2)
104 : fswthrun_idf, & ! nir dif sw radiation through ice bot (W/m**2)
105 : Urefn ! air speed reference level (m/s)
106 :
107 : ! cumulative fluxes
108 : real (kind=dbl_kind), optional, intent(inout) :: &
109 : strairxT, & ! air/ice zonal strss, (N/m**2)
110 : strairyT, & ! air/ice merdnl strss, (N/m**2)
111 : Cdn_atm_ratio, & ! ratio of total drag over neutral drag
112 : fsurf , & ! net heat flux to top surface (W/m**2)
113 : fcondtop, & ! downward cond flux at top sfc (W/m**2)
114 : fcondbot, & ! downward cond flux at bottom sfc (W/m**2)
115 : fsens , & ! sensible heat flx (W/m**2)
116 : flat , & ! latent heat flx (W/m**2)
117 : fswabs , & ! shortwave absorbed heat flx (W/m**2)
118 : flwout , & ! upwd lw emitted heat flx (W/m**2)
119 : evap , & ! evaporation (kg/m2/s)
120 : evaps , & ! evaporation over snow (kg/m2/s)
121 : evapi , & ! evaporation over ice (kg/m2/s)
122 : Tref , & ! air tmp reference level (K)
123 : Qref , & ! air sp hum reference level (kg/kg)
124 : fresh , & ! fresh water flux to ocean (kg/m2/s)
125 : fsalt , & ! salt flux to ocean (kg/m2/s)
126 : fhocn , & ! actual ocn/ice heat flx (W/m**2)
127 : fswthru , & ! sw radiation through ice bot (W/m**2)
128 : meltt , & ! top ice melt (m)
129 : meltb , & ! bottom ice melt (m)
130 : melts , & ! snow melt (m)
131 : meltsliq, & ! mass of snow melt (kg/m^2)
132 : congel , & ! congelation ice growth (m)
133 : snoice , & ! snow-ice growth (m)
134 : fswthru_vdr, & ! vis dir sw radiation through ice bot (W/m**2)
135 : fswthru_vdf, & ! vis dif sw radiation through ice bot (W/m**2)
136 : fswthru_idr, & ! nir dir sw radiation through ice bot (W/m**2)
137 : fswthru_idf, & ! nir dif sw radiation through ice bot (W/m**2)
138 : dsnow, & ! change in snow depth (m)
139 : Uref ! air speed reference level (m/s)
140 :
141 : real (kind=dbl_kind), dimension(:), intent(in), optional :: &
142 : Qrefn_iso, & ! isotope air sp hum ref level (kg/kg)
143 : fiso_ocnn, & ! isotope fluxes to ocean (kg/m2/s)
144 : fiso_evapn ! isotope evaporation (kg/m2/s)
145 :
146 : real (kind=dbl_kind), dimension(:), intent(inout), optional :: &
147 : Qref_iso, & ! isotope air sp hum ref level (kg/kg)
148 : fiso_ocn, & ! isotope fluxes to ocean (kg/m2/s)
149 : fiso_evap ! isotope evaporation (kg/m2/s)
150 :
151 : character(len=*),parameter :: subname='(merge_fluxes)'
152 :
153 : !-----------------------------------------------------------------
154 : ! Merge fluxes
155 : ! NOTE: The albedo is aggregated only in cells where ice exists
156 : ! and (for the delta-Eddington scheme) where the sun is above
157 : ! the horizon.
158 : !-----------------------------------------------------------------
159 :
160 : ! atmo fluxes
161 :
162 8777938 : if (present(strairxn) .and. present(strairxT)) &
163 8777938 : strairxT = strairxT + strairxn * aicen
164 8777938 : if (present(strairyn) .and. present(strairyT)) &
165 8777938 : strairyT = strairyT + strairyn * aicen
166 8777938 : if (present(Cdn_atm_ratio_n) .and. present(Cdn_atm_ratio)) &
167 : Cdn_atm_ratio = Cdn_atm_ratio + &
168 8777938 : Cdn_atm_ratio_n * aicen
169 8777938 : if (present(fsurfn) .and. present(fsurf)) &
170 8777938 : fsurf = fsurf + fsurfn * aicen
171 8777938 : if (present(fcondtopn) .and. present(fcondtop)) &
172 8777938 : fcondtop = fcondtop + fcondtopn * aicen
173 8777938 : if (present(fcondbotn) .and. present(fcondbot)) &
174 8777938 : fcondbot = fcondbot + fcondbotn * aicen
175 8777938 : if (present(fsensn) .and. present(fsens)) &
176 8777938 : fsens = fsens + fsensn * aicen
177 8777938 : if (present(flatn) .and. present(flat)) &
178 8777938 : flat = flat + flatn * aicen
179 8777938 : if (present(fswabsn) .and. present(fswabs)) &
180 8777938 : fswabs = fswabs + fswabsn * aicen
181 8777938 : if (present(flwoutn) .and. present(flwout) .and. present(flw)) &
182 : flwout = flwout &
183 8777938 : + (flwoutn - (c1-emissivity)*flw) * aicen
184 8777938 : if (present(evapn) .and. present(evap)) &
185 8777938 : evap = evap + evapn * aicen
186 8777938 : if (present(evapsn) .and. present(evaps)) &
187 8777938 : evaps = evaps + evapsn * aicen
188 8777938 : if (present(evapin) .and. present(evapi)) &
189 8777938 : evapi = evapi + evapin * aicen
190 8777938 : if (present(Trefn) .and. present(Tref)) &
191 8777938 : Tref = Tref + Trefn * aicen
192 8777938 : if (present(Qrefn) .and. present(Qref)) &
193 8777938 : Qref = Qref + Qrefn * aicen
194 :
195 : ! Isotopes
196 8777938 : if (tr_iso) then
197 225452 : if (present(Qrefn_iso) .and. present(Qref_iso)) then
198 901808 : Qref_iso (:) = Qref_iso (:) + Qrefn_iso (:) * aicen
199 : endif
200 225452 : if (present(fiso_ocnn) .and. present(fiso_ocn)) then
201 901808 : fiso_ocn (:) = fiso_ocn (:) + fiso_ocnn (:) * aicen
202 : endif
203 225452 : if (present(fiso_evapn) .and. present(fiso_evap)) then
204 901808 : fiso_evap(:) = fiso_evap(:) + fiso_evapn(:) * aicen
205 : endif
206 : endif
207 :
208 : ! ocean fluxes
209 8777938 : if (present(Urefn) .and. present(Uref)) then
210 8777938 : Uref = Uref + Urefn * aicen
211 : endif
212 :
213 8777938 : if (present(freshn) .and. present(fresh)) &
214 8777938 : fresh = fresh + freshn * aicen
215 8777938 : if (present(fsaltn) .and. present(fsalt)) &
216 8777938 : fsalt = fsalt + fsaltn * aicen
217 8777938 : if (present(fhocnn) .and. present(fhocn)) &
218 8777938 : fhocn = fhocn + fhocnn * aicen
219 8777938 : if (present(fswthrun) .and. present(fswthru)) &
220 8777938 : fswthru = fswthru + fswthrun * aicen
221 :
222 8777938 : if (present(fswthrun_vdr) .and. present(fswthru_vdr)) &
223 8777938 : fswthru_vdr = fswthru_vdr + fswthrun_vdr * aicen
224 8777938 : if (present(fswthrun_vdf) .and. present(fswthru_vdf)) &
225 8777938 : fswthru_vdf = fswthru_vdf + fswthrun_vdf * aicen
226 8777938 : if (present(fswthrun_idr) .and. present(fswthru_idr)) &
227 8777938 : fswthru_idr = fswthru_idr + fswthrun_idr * aicen
228 8777938 : if (present(fswthrun_idf) .and. present(fswthru_idf)) &
229 8777938 : fswthru_idf = fswthru_idf + fswthrun_idf * aicen
230 :
231 : ! ice/snow thickness
232 :
233 8777938 : if (present(melttn) .and. present(meltt)) &
234 8777938 : meltt = meltt + melttn * aicen
235 8777938 : if (present(meltbn) .and. present(meltb)) &
236 8777938 : meltb = meltb + meltbn * aicen
237 8777938 : if (present(meltsn) .and. present(melts)) &
238 8777938 : melts = melts + meltsn * aicen
239 8777938 : if (snwgrain) then
240 885289 : if (present(meltsliqn) .and. present(meltsliq)) &
241 885289 : meltsliq = meltsliq + meltsliqn * aicen
242 : endif
243 8777938 : if (present(dsnown) .and. present(dsnow)) then
244 8777938 : dsnow = dsnow + dsnown * aicen
245 : endif
246 8777938 : if (present(congeln) .and. present(congel)) &
247 8777938 : congel = congel + congeln * aicen
248 8777938 : if (present(snoicen) .and. present(snoice)) &
249 8777938 : snoice = snoice + snoicen * aicen
250 :
251 8777938 : end subroutine merge_fluxes
252 :
253 : !=======================================================================
254 :
255 : ! If model is not calculating surface temperature, set the surface
256 : ! flux values using values read in from forcing data or supplied via
257 : ! coupling (stored in ice_flux).
258 : !
259 : ! If CICE is running in NEMO environment, convert fluxes from GBM values
260 : ! to per unit ice area values. If model is not running in NEMO environment,
261 : ! the forcing is supplied as per unit ice area values.
262 : !
263 : ! authors Alison McLaren, Met Office
264 :
265 0 : subroutine set_sfcflux (aicen, &
266 : flatn_f, &
267 : fsensn_f, &
268 : fsurfn_f, &
269 : fcondtopn_f, &
270 : flatn, &
271 : fsensn, &
272 : fsurfn, &
273 : fcondtopn)
274 :
275 : ! ice state variables
276 : real (kind=dbl_kind), intent(in) :: &
277 : aicen , & ! concentration of ice
278 : flatn_f , & ! latent heat flux (W/m^2)
279 : fsensn_f , & ! sensible heat flux (W/m^2)
280 : fsurfn_f , & ! net flux to top surface, not including fcondtopn
281 : fcondtopn_f ! downward cond flux at top surface (W m-2)
282 :
283 : real (kind=dbl_kind), intent(out):: &
284 : flatn , & ! latent heat flux (W/m^2)
285 : fsensn , & ! sensible heat flux (W/m^2)
286 : fsurfn , & ! net flux to top surface, not including fcondtopn
287 : fcondtopn ! downward cond flux at top surface (W m-2)
288 :
289 : ! local variables
290 :
291 : real (kind=dbl_kind) :: &
292 : raicen ! 1 or 1/aicen
293 :
294 : logical (kind=log_kind) :: &
295 : extreme_flag ! flag for extreme forcing values
296 :
297 : logical (kind=log_kind), parameter :: &
298 : extreme_test=.false. ! test and write out extreme forcing data
299 :
300 : character(len=*),parameter :: subname='(set_sfcflux)'
301 :
302 0 : raicen = c1
303 :
304 : #ifdef CICE_IN_NEMO
305 : !----------------------------------------------------------------------
306 : ! Convert fluxes from GBM values to per ice area values when
307 : ! running in NEMO environment. (When in standalone mode, fluxes
308 : ! are input as per ice area.)
309 : !----------------------------------------------------------------------
310 : raicen = c1 / aicen
311 : #endif
312 0 : fsurfn = fsurfn_f*raicen
313 0 : fcondtopn= fcondtopn_f*raicen
314 0 : flatn = flatn_f*raicen
315 0 : fsensn = fsensn_f*raicen
316 :
317 : !----------------------------------------------------------------
318 : ! Flag up any extreme fluxes
319 : !---------------------------------------------------------------
320 :
321 : if (extreme_test) then
322 : extreme_flag = .false.
323 :
324 : if (fcondtopn < -100.0_dbl_kind &
325 : .or. fcondtopn > 20.0_dbl_kind) then
326 : extreme_flag = .true.
327 : endif
328 :
329 : if (fsurfn < -100.0_dbl_kind &
330 : .or. fsurfn > 80.0_dbl_kind) then
331 : extreme_flag = .true.
332 : endif
333 :
334 : if (flatn < -20.0_dbl_kind &
335 : .or. flatn > 20.0_dbl_kind) then
336 : extreme_flag = .true.
337 : endif
338 :
339 : if (extreme_flag) then
340 :
341 : if (fcondtopn < -100.0_dbl_kind &
342 : .or. fcondtopn > 20.0_dbl_kind) then
343 : write(warnstr,*) subname, &
344 : 'Extreme forcing: -100 > fcondtopn > 20'
345 : call icepack_warnings_add(warnstr)
346 : write(warnstr,*) subname, &
347 : 'aicen,fcondtopn = ', &
348 : aicen,fcondtopn
349 : call icepack_warnings_add(warnstr)
350 : endif
351 :
352 : if (fsurfn < -100.0_dbl_kind &
353 : .or. fsurfn > 80.0_dbl_kind) then
354 : write(warnstr,*) subname, &
355 : 'Extreme forcing: -100 > fsurfn > 40'
356 : call icepack_warnings_add(warnstr)
357 : write(warnstr,*) subname, &
358 : 'aicen,fsurfn = ', &
359 : aicen,fsurfn
360 : call icepack_warnings_add(warnstr)
361 : endif
362 :
363 : if (flatn < -20.0_dbl_kind &
364 : .or. flatn > 20.0_dbl_kind) then
365 : write(warnstr,*) subname, &
366 : 'Extreme forcing: -20 > flatn > 20'
367 : call icepack_warnings_add(warnstr)
368 : write(warnstr,*) subname, &
369 : 'aicen,flatn = ', &
370 : aicen,flatn
371 : call icepack_warnings_add(warnstr)
372 : endif
373 :
374 : endif ! extreme_flag
375 : endif ! extreme_test
376 :
377 0 : end subroutine set_sfcflux
378 :
379 : !=======================================================================
380 :
381 : end module icepack_flux
382 :
383 : !=======================================================================
|