Line data Source code
1 : !=======================================================================
2 :
3 : ! Ocean boundary interface
4 :
5 : module icepack_ocean
6 :
7 : use icepack_kinds
8 : use icepack_parameters, only: c0, c1, c1000
9 : use icepack_parameters, only: Tffresh, stefan_boltzmann, Lvap, cprho
10 : use icepack_warnings, only: warnstr, icepack_warnings_add
11 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
12 :
13 : implicit none
14 :
15 : private
16 : public :: icepack_ocn_mixed_layer
17 :
18 : !=======================================================================
19 :
20 : contains
21 :
22 : !=======================================================================
23 : !=======================================================================
24 : !autodocument_start icepack_ocn_mixed_layer
25 : ! Compute the mixed layer heat balance and update the SST.
26 : ! Compute the energy available to freeze or melt ice.
27 : ! NOTE: SST changes due to fluxes through the ice are computed in
28 : ! icepack_therm_vertical.
29 :
30 945373128 : subroutine icepack_ocn_mixed_layer (alvdr_ocn, swvdr, &
31 : alidr_ocn, swidr, &
32 : alvdf_ocn, swvdf, &
33 : alidf_ocn, swidf, &
34 : sst, flwout_ocn, &
35 : fsens_ocn, shcoef, &
36 : flat_ocn, lhcoef, &
37 : evap_ocn, flw, &
38 : delt, delq, &
39 : aice, fhocn, &
40 : fswthru, hmix, &
41 : Tf, qdp, &
42 : frzmlt, dt)
43 :
44 : real (kind=dbl_kind), intent(in) :: &
45 : alvdr_ocn , & ! visible, direct (fraction)
46 : alidr_ocn , & ! near-ir, direct (fraction)
47 : alvdf_ocn , & ! visible, diffuse (fraction)
48 : alidf_ocn , & ! near-ir, diffuse (fraction)
49 : swvdr , & ! sw down, visible, direct (W/m^2)
50 : swvdf , & ! sw down, visible, diffuse (W/m^2)
51 : swidr , & ! sw down, near IR, direct (W/m^2)
52 : swidf , & ! sw down, near IR, diffuse (W/m^2)
53 : flw , & ! incoming longwave radiation (W/m^2)
54 : Tf , & ! freezing temperature (C)
55 : hmix , & ! mixed layer depth (m)
56 : delt , & ! potential temperature difference (K)
57 : delq , & ! specific humidity difference (kg/kg)
58 : shcoef , & ! transfer coefficient for sensible heat
59 : lhcoef , & ! transfer coefficient for latent heat
60 : fhocn , & ! net heat flux to ocean (W/m^2)
61 : fswthru , & ! shortwave penetrating to ocean (W/m^2)
62 : aice , & ! ice area fraction
63 : dt ! time step (s)
64 :
65 : real (kind=dbl_kind), intent(inout) :: &
66 : flwout_ocn, & ! outgoing longwave radiation (W/m^2)
67 : fsens_ocn , & ! sensible heat flux (W/m^2)
68 : flat_ocn , & ! latent heat flux (W/m^2)
69 : evap_ocn , & ! evaporative water flux (kg/m^2/s)
70 : qdp , & ! deep ocean heat flux (W/m^2), negative upward
71 : sst , & ! sea surface temperature (C)
72 : frzmlt ! freezing/melting potential (W/m^2)
73 :
74 : !autodocument_end
75 :
76 : ! local variables
77 :
78 : real (kind=dbl_kind), parameter :: &
79 : frzmlt_max = c1000 ! max magnitude of frzmlt (W/m^2)
80 :
81 : real (kind=dbl_kind) :: &
82 59897256 : TsfK , & ! surface temperature (K)
83 59897256 : swabs ! surface absorbed shortwave heat flux (W/m^2)
84 :
85 : character(len=*),parameter :: subname='(icepack_ocn_mixed_layer)'
86 :
87 : ! shortwave radiative flux
88 : swabs = (c1-alvdr_ocn) * swvdr + (c1-alidr_ocn) * swidr &
89 945373128 : + (c1-alvdf_ocn) * swvdf + (c1-alidf_ocn) * swidf
90 :
91 : ! ocean surface temperature in Kelvin
92 945373128 : TsfK = sst + Tffresh
93 :
94 : ! longwave radiative flux
95 945373128 : flwout_ocn = -stefan_boltzmann * TsfK**4
96 :
97 : ! downward latent and sensible heat fluxes
98 945373128 : fsens_ocn = shcoef * delt
99 945373128 : flat_ocn = lhcoef * delq
100 945373128 : evap_ocn = -flat_ocn / Lvap
101 :
102 : ! Compute sst change due to exchange with atm/ice above
103 : sst = sst + dt * ( &
104 : (fsens_ocn + flat_ocn + flwout_ocn + flw + swabs) * (c1-aice) &
105 : + fhocn + fswthru) & ! these are *aice
106 945373128 : / (cprho*hmix)
107 :
108 : ! adjust qdp if cooling of mixed layer would occur when sst <= Tf
109 945373128 : if (sst <= Tf .and. qdp > c0) qdp = c0
110 :
111 : ! computed T change due to exchange with deep layers:
112 945373128 : sst = sst - qdp*dt/(cprho*hmix)
113 :
114 : ! compute potential to freeze or melt ice
115 945373128 : frzmlt = (Tf-sst)*cprho*hmix/dt
116 945373128 : frzmlt = min(max(frzmlt,-frzmlt_max),frzmlt_max)
117 :
118 : ! if sst is below freezing, reset sst to Tf
119 945373128 : if (sst <= Tf) sst = Tf
120 :
121 945373128 : end subroutine icepack_ocn_mixed_layer
122 :
123 : !=======================================================================
124 :
125 : end module icepack_ocean
126 :
127 : !=======================================================================
|