-
Notifications
You must be signed in to change notification settings - Fork 138
/
Copy pathicepack_therm_shared.F90
581 lines (418 loc) · 18.4 KB
/
icepack_therm_shared.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
!=========================================================================
!
! Shared thermo variables, subroutines
!
! authors: Elizabeth C. Hunke, LANL
module icepack_therm_shared
use icepack_kinds
use icepack_parameters, only: c0, c1, c2, c4, p5, pi, puny, Tocnfrz
use icepack_parameters, only: cp_ocn, cp_ice, rhoi, rhos, Tffresh, TTTice, qqqice
use icepack_parameters, only: stefan_boltzmann, emissivity, Lfresh, Tsmelt
use icepack_parameters, only: saltmax, min_salin, depressT
use icepack_parameters, only: ktherm, tfrz_option
use icepack_parameters, only: calc_Tsfc
use icepack_warnings, only: warnstr, icepack_warnings_add
use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
use icepack_mushy_physics, only: icepack_enthalpy_mush
use icepack_mushy_physics, only: temperature_snow
use icepack_mushy_physics, only: icepack_mushy_temperature_mush
use icepack_mushy_physics, only: liquidus_temperature_mush
implicit none
private
public :: calculate_Tin_from_qin, &
surface_heat_flux, &
dsurface_heat_flux_dTsf, &
icepack_init_thermo, &
icepack_salinity_profile, &
icepack_init_trcr, &
icepack_ice_temperature, &
icepack_snow_temperature, &
icepack_liquidus_temperature, &
icepack_sea_freezing_temperature, &
adjust_enthalpy
real (kind=dbl_kind), parameter, public :: &
ferrmax = 1.0e-3_dbl_kind ! max allowed energy flux error (W m-2)
! recommend ferrmax < 0.01 W m-2
real (kind=dbl_kind), parameter, public :: &
Tmin = -100.0_dbl_kind ! min allowed internal temperature (deg C)
logical (kind=log_kind), public :: &
l_brine ! if true, treat brine pocket effects
!=======================================================================
contains
!=======================================================================
!
! Compute the internal ice temperatures from enthalpy using
! quadratic formula
function calculate_Tin_from_qin (qin, Tmltk) &
result(Tin)
real (kind=dbl_kind), intent(in) :: &
qin , & ! enthalpy
Tmltk ! melting temperature at one level
real (kind=dbl_kind) :: &
Tin ! internal temperature
! local variables
real (kind=dbl_kind) :: &
aa1,bb1,cc1 ! quadratic solvers
character(len=*),parameter :: subname='(calculate_Tin_from_qin)'
if (l_brine) then
aa1 = cp_ice
bb1 = (cp_ocn-cp_ice)*Tmltk - qin/rhoi - Lfresh
cc1 = Lfresh * Tmltk
Tin = min((-bb1 - sqrt(bb1*bb1 - c4*aa1*cc1)) / &
(c2*aa1),Tmltk)
else ! fresh ice
Tin = (Lfresh + qin/rhoi) / cp_ice
endif
end function calculate_Tin_from_qin
!=======================================================================
! Surface heat flux
!=======================================================================
! heat flux into ice
subroutine surface_heat_flux(Tsf, fswsfc, &
rhoa, flw, &
potT, Qa, &
shcoef, lhcoef, &
flwoutn, fsensn, &
flatn, fsurfn)
! input surface temperature
real(kind=dbl_kind), intent(in) :: &
Tsf ! ice/snow surface temperature (C)
! input variables
real(kind=dbl_kind), intent(in) :: &
fswsfc , & ! SW absorbed at ice/snow surface (W m-2)
rhoa , & ! air density (kg/m^3)
flw , & ! incoming longwave radiation (W/m^2)
potT , & ! air potential temperature (K)
Qa , & ! specific humidity (kg/kg)
shcoef , & ! transfer coefficient for sensible heat
lhcoef ! transfer coefficient for latent heat
! output
real(kind=dbl_kind), intent(out) :: &
fsensn , & ! surface downward sensible heat (W m-2)
flatn , & ! surface downward latent heat (W m-2)
flwoutn , & ! upward LW at surface (W m-2)
fsurfn ! net flux to top surface, excluding fcondtopn
! local variables
real(kind=dbl_kind) :: &
TsfK , & ! ice/snow surface temperature (K)
Qsfc , & ! saturated surface specific humidity (kg/kg)
qsat , & ! the saturation humidity of air (kg/m^3)
flwdabs , & ! downward longwave absorbed heat flx (W/m^2)
tmpvar ! 1/TsfK
character(len=*),parameter :: subname='(surface_heat_flux)'
! ice surface temperature in Kelvin
TsfK = Tsf + Tffresh
! TsfK = max(Tsf + Tffresh, c1)
tmpvar = c1/TsfK
! saturation humidity
qsat = qqqice * exp(-TTTice*tmpvar)
Qsfc = qsat / rhoa
! longwave radiative flux
flwdabs = emissivity * flw
flwoutn = -emissivity * stefan_boltzmann * TsfK**4
! downward latent and sensible heat fluxes
fsensn = shcoef * (potT - TsfK)
flatn = lhcoef * (Qa - Qsfc)
! combine fluxes
fsurfn = fswsfc + flwdabs + flwoutn + fsensn + flatn
end subroutine surface_heat_flux
!=======================================================================
subroutine dsurface_heat_flux_dTsf(Tsf, rhoa, &
shcoef, lhcoef, &
dfsurfn_dTsf, dflwoutn_dTsf, &
dfsensn_dTsf, dflatn_dTsf)
! input surface temperature
real(kind=dbl_kind), intent(in) :: &
Tsf ! ice/snow surface temperature (C)
! input variables
real(kind=dbl_kind), intent(in) :: &
rhoa , & ! air density (kg/m^3)
shcoef , & ! transfer coefficient for sensible heat
lhcoef ! transfer coefficient for latent heat
! output
real(kind=dbl_kind), intent(out) :: &
dfsurfn_dTsf ! derivative of net flux to top surface, excluding fcondtopn
real(kind=dbl_kind), intent(out) :: &
dflwoutn_dTsf , & ! derivative of longwave flux wrt surface temperature
dfsensn_dTsf , & ! derivative of sensible heat flux wrt surface temperature
dflatn_dTsf ! derivative of latent heat flux wrt surface temperature
! local variables
real(kind=dbl_kind) :: &
TsfK , & ! ice/snow surface temperature (K)
dQsfc_dTsf , & ! saturated surface specific humidity (kg/kg)
qsat , & ! the saturation humidity of air (kg/m^3)
tmpvar ! 1/TsfK
character(len=*),parameter :: subname='(dsurface_heat_flux_dTsf)'
! ice surface temperature in Kelvin
! TsfK = max(Tsf + Tffresh, c1)
TsfK = Tsf + Tffresh
tmpvar = c1/TsfK
! saturation humidity
qsat = qqqice * exp(-TTTice*tmpvar)
dQsfc_dTsf = TTTice * tmpvar * tmpvar * (qsat / rhoa)
! longwave radiative flux
dflwoutn_dTsf = -emissivity * stefan_boltzmann * c4*TsfK**3
! downward latent and sensible heat fluxes
dfsensn_dTsf = -shcoef
dflatn_dTsf = -lhcoef * dQsfc_dTsf
! combine fluxes
dfsurfn_dTsf = dflwoutn_dTsf + dfsensn_dTsf + dflatn_dTsf
end subroutine dsurface_heat_flux_dTsf
!=======================================================================
!autodocument_start icepack_init_thermo
! Initialize the vertical profile of ice salinity and melting temperature.
!
! authors: C. M. Bitz, UW
! William H. Lipscomb, LANL
subroutine icepack_init_thermo(nilyr, sprofile)
integer (kind=int_kind), intent(in) :: &
nilyr ! number of ice layers
real (kind=dbl_kind), dimension(:), intent(out) :: &
sprofile ! vertical salinity profile
!autodocument_end
integer (kind=int_kind) :: k ! ice layer index
real (kind=dbl_kind) :: zn ! normalized ice thickness
character(len=*),parameter :: subname='(icepack_init_thermo)'
!-----------------------------------------------------------------
! Determine l_brine based on saltmax.
! Thermodynamic solver will not converge if l_brine is true and
! saltmax is close to zero.
! Set l_brine to false for zero layer thermodynamics
!-----------------------------------------------------------------
l_brine = .false.
if (saltmax > min_salin) l_brine = .true.
!-----------------------------------------------------------------
! Prescibe vertical profile of salinity and melting temperature.
! Note this profile is only used for BL99 thermodynamics.
!-----------------------------------------------------------------
if (l_brine) then
do k = 1, nilyr
zn = (real(k,kind=dbl_kind)-p5) / &
real(nilyr,kind=dbl_kind)
! sprofile(k)=(saltmax/c2)*(c1-cos(pi*zn**(nsal/(msal+zn))))
sprofile(k)=icepack_salinity_profile(zn)
sprofile(k) = max(sprofile(k), min_salin)
enddo ! k
sprofile(nilyr+1) = saltmax
else ! .not. l_brine
do k = 1, nilyr+1
sprofile(k) = c0
enddo
endif ! l_brine
end subroutine icepack_init_thermo
!=======================================================================
!autodocument_start icepack_salinity_profile
! Initial salinity profile
!
! authors: C. M. Bitz, UW
! William H. Lipscomb, LANL
function icepack_salinity_profile(zn) result(salinity)
real(kind=dbl_kind), intent(in) :: &
zn ! depth
real(kind=dbl_kind) :: &
salinity ! initial salinity profile
!autodocument_end
real (kind=dbl_kind), parameter :: &
nsal = 0.407_dbl_kind, &
msal = 0.573_dbl_kind
character(len=*),parameter :: subname='(icepack_init_thermo)'
salinity = (saltmax/c2)*(c1-cos(pi*zn**(nsal/(msal+zn))))
end function icepack_salinity_profile
!=======================================================================
!autodocument_start icepack_init_trcr
!
subroutine icepack_init_trcr(Tair, Tf, &
Sprofile, Tprofile, &
Tsfc, &
nilyr, nslyr, &
qin, qsn)
integer (kind=int_kind), intent(in) :: &
nilyr, & ! number of ice layers
nslyr ! number of snow layers
real (kind=dbl_kind), intent(in) :: &
Tair, & ! air temperature (K)
Tf ! freezing temperature (C)
real (kind=dbl_kind), dimension(:), intent(in) :: &
Sprofile, & ! vertical salinity profile (ppt)
Tprofile ! vertical temperature profile (C)
real (kind=dbl_kind), intent(out) :: &
Tsfc ! surface temperature (C)
real (kind=dbl_kind), dimension(:), intent(out) :: &
qin, & ! ice enthalpy profile (J/m3)
qsn ! snow enthalpy profile (J/m3)
!autodocument_end
! local variables
integer (kind=int_kind) :: k
real (kind=dbl_kind) :: &
slope, Ti
character(len=*),parameter :: subname='(icepack_init_trcr)'
! surface temperature
Tsfc = Tf ! default
if (calc_Tsfc) Tsfc = min(Tsmelt, Tair - Tffresh) ! deg C
! ice enthalpy
do k = 1, nilyr
! assume linear temp profile and compute enthalpy
slope = Tf - Tsfc
Ti = Tsfc + slope*(real(k,kind=dbl_kind)-p5) &
/real(nilyr,kind=dbl_kind)
if (ktherm == 2) then
qin(k) = icepack_enthalpy_mush(Ti, Sprofile(k))
else
qin(k) = -(rhoi * (cp_ice*(Tprofile(k)-Ti) &
+ Lfresh*(c1-Tprofile(k)/Ti) - cp_ocn*Tprofile(k)))
endif
enddo ! nilyr
! snow enthalpy
do k = 1, nslyr
Ti = min(c0, Tsfc)
qsn(k) = -rhos*(Lfresh - cp_ice*Ti)
enddo ! nslyr
end subroutine icepack_init_trcr
!=======================================================================
!autodocument_start icepack_liquidus_temperature
! compute liquidus temperature
function icepack_liquidus_temperature(Sin) result(Tmlt)
real(dbl_kind), intent(in) :: Sin
real(dbl_kind) :: Tmlt
!autodocument_end
character(len=*),parameter :: subname='(icepack_liquidus_temperature)'
if (ktherm == 2) then
Tmlt = liquidus_temperature_mush(Sin)
else
Tmlt = -depressT * Sin
endif
end function icepack_liquidus_temperature
!=======================================================================
!autodocument_start icepack_sea_freezing_temperature
! compute ocean freezing temperature
function icepack_sea_freezing_temperature(sss) result(Tf)
real(dbl_kind), intent(in) :: sss
real(dbl_kind) :: Tf
!autodocument_end
character(len=*),parameter :: subname='(icepack_sea_freezing_temperature)'
if (trim(tfrz_option) == 'mushy') then
Tf = icepack_liquidus_temperature(sss) ! deg C
elseif (trim(tfrz_option) == 'linear_salt') then
Tf = -depressT * sss ! deg C
elseif (trim(tfrz_option) == 'constant') then
Tf = Tocnfrz
elseif (trim(tfrz_option) == 'minus1p8') then
Tf = -1.8_dbl_kind
else
call icepack_warnings_add(subname//' tfrz_option unsupported: '//trim(tfrz_option))
call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
endif
end function icepack_sea_freezing_temperature
!=======================================================================
!autodocument_start icepack_ice_temperature
! compute ice temperature
function icepack_ice_temperature(qin, Sin) result(Tin)
real(kind=dbl_kind), intent(in) :: qin, Sin
real(kind=dbl_kind) :: Tin
!autodocument_end
real(kind=dbl_kind) :: Tmlts
character(len=*),parameter :: subname='(icepack_ice_temperature)'
if (ktherm == 2) then
Tin = icepack_mushy_temperature_mush(qin, Sin)
else
Tmlts = -depressT * Sin
Tin = calculate_Tin_from_qin(qin,Tmlts)
endif
end function icepack_ice_temperature
!=======================================================================
!autodocument_start icepack_snow_temperature
! compute snow temperature
function icepack_snow_temperature(qin) result(Tsn)
real(kind=dbl_kind), intent(in) :: qin
real(kind=dbl_kind) :: Tsn
!autodocument_end
character(len=*),parameter :: subname='(icepack_snow_temperature)'
if (ktherm == 2) then
Tsn = temperature_snow(qin)
else
Tsn = (Lfresh + qin/rhos)/cp_ice
endif
end function icepack_snow_temperature
!=======================================================================
!
! Conserving energy, compute the new enthalpy of equal-thickness ice
! or snow layers.
!
! authors William H. Lipscomb, LANL
! C. M. Bitz, UW
subroutine adjust_enthalpy (nlyr, &
z1, z2, &
hlyr, hn, &
qn)
integer (kind=int_kind), intent(in) :: &
nlyr ! number of layers (nilyr or nslyr)
real (kind=dbl_kind), dimension (:), intent(in) :: &
z1 , & ! interface depth for old, unequal layers (m)
z2 ! interface depth for new, equal layers (m)
real (kind=dbl_kind), intent(in) :: &
hlyr ! new layer thickness (m)
real (kind=dbl_kind), intent(in) :: &
hn ! total thickness (m)
real (kind=dbl_kind), dimension (:), intent(inout) :: &
qn ! layer quantity (enthalpy, salinity...)
! local variables
integer (kind=int_kind) :: &
k, k1, k2 ! vertical indices
real (kind=dbl_kind) :: &
hovlp ! overlap between old and new layers (m)
real (kind=dbl_kind) :: &
rhlyr, & ! 1./hlyr
qtot ! total h*q in the column
real (kind=dbl_kind), dimension (nlyr) :: &
hq ! h * q for a layer
character(len=*),parameter :: subname='(adjust_enthalpy)'
!-----------------------------------------------------------------
! Compute reciprocal layer thickness.
!-----------------------------------------------------------------
rhlyr = c0
if (hn > puny) then
rhlyr = c1 / hlyr
!-----------------------------------------------------------------
! Compute h*q for new layers (k2) given overlap with old layers (k1)
!-----------------------------------------------------------------
do k2 = 1, nlyr
hq(k2) = c0
enddo ! k
k1 = 1
k2 = 1
do while (k1 <= nlyr .and. k2 <= nlyr)
hovlp = min (z1(k1+1), z2(k2+1)) &
- max (z1(k1), z2(k2))
hovlp = max (hovlp, c0)
hq(k2) = hq(k2) + hovlp*qn(k1)
if (z1(k1+1) > z2(k2+1)) then
k2 = k2 + 1
else
k1 = k1 + 1
endif
enddo ! while
!-----------------------------------------------------------------
! Compute new enthalpies.
!-----------------------------------------------------------------
do k = 1, nlyr
qn(k) = hq(k) * rhlyr
enddo ! k
else
qtot = c0
do k = 1, nlyr
qtot = qtot + qn(k) * (z1(k+1)-z1(k))
enddo
if (hn > c0) then
do k = 1, nlyr
qn(k) = qtot/hn
enddo
else
do k = 1, nlyr
qn(k) = c0
enddo
endif
endif
end subroutine adjust_enthalpy
!=======================================================================
end module icepack_therm_shared
!=======================================================================