My Project
ice_ocean.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! CVS VERSION INFORMATION
14 ! $Id$
15 ! $Name$
16 ! $Revision$
17 !/===========================================================================/
18 
19 !=======================================================================
20 !BOP
21 !
22 ! !MODULE: ice_ocean - ocean mixed layer internal to sea ice model
23 !
24 ! !DESCRIPTION:
25 !
26 ! Ocean mixed layer calculation (internal to sea ice model).
27 ! Allows heat storage in ocean for uncoupled runs.
28 !
29 ! !REVISION HISTORY:
30 !
31 ! authors: John Weatherly, CRREL
32 ! C.M. Bitz, UW
33 ! Elizabeth C. Hunke, LANL
34 ! Bruce P. Briegleb, NCAR
35 ! William H. Lipscomb, LANL
36 !
37 ! !INTERFACE:
38 !
39  module ice_ocean
40 !
41 ! !USES:
42 !
43  use ice_kinds_mod
44  use ice_constants
45 ! use ice_calendar, only: dt
46  use ice_calendar, only: dtice
47 !
48 !EOP
49 !
50  implicit none
51  save
52 
53  logical (kind=log_kind) :: &
54  oceanmixed_ice ! if true, use ocean mixed layer
55 
56 !=======================================================================
57 
58  contains
59 
60 !=======================================================================
61 !BOP
62 !
63 ! !ROUTINE: mixed_layer - compute SST and freeze/melt potential
64 !
65 ! !DESCRIPTION:
66 !
67 ! Compute the mixed layer heat balance and update the SST.
68 ! Compute the energy available to freeze or melt ice.
69 ! NOTE: SST changes due to fluxes through the ice are computed in
70 ! ice_therm_vertical.
71 !
72 ! !REVISION HISTORY:
73 !
74 ! authors: John Weatherly, CRREL
75 ! C.M. Bitz, UW
76 ! Elizabeth C. Hunke, LANL
77 ! Bruce P. Briegleb, NCAR
78 ! William H. Lipscomb, LANL
79 !
80 ! !INTERFACE:
81 !
82  subroutine mixed_layer
83 !
84 ! !USES:
85 !
86  use ice_flux
87  use ice_grid, only: tmask
88  use ice_atmo
89  use ice_state
90  use ice_albedo
91 !
92 ! !INPUT/OUTPUT PARAMETERS:
93 !
94 !EOP
95 !
96  integer (kind=int_kind) :: &
97  i, j ! horizontal indices
98 
99  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi) ::&
100  delt & ! potential temperature difference (K)
101  , delq & ! specific humidity difference (kg/kg)
102  , dummy1, dummy2, dummy3, dummy4 ! dummy arrays
103 
104  real (kind=dbl_kind) :: &
105  tsfk &! surface temperature (K)
106  , fsh &! sensible heat flux (W/m^2)
107  , flh &! latent heat flux (W/m^2)
108  , swabs &! surface absorbed shortwave heat flux (W/m^2)
109  , flwup &! long-wave upward heat flux (W/m^2)
110  , ft &! fraction reduction of positive qdp
111  , dtcprho ! dt/cp_ocn*rhow
112 
113  dtcprho = dtice/(cp_ocn*rhow)
114 
115  call atmo_boundary_layer (1, 'ocn', sst, &
116  dummy1, dummy2, dummy3, dummy4, delt, delq)
117 
118  do j = jlo,jhi
119  do i = ilo,ihi
120  if (tmask(i,j)) then
121 ! if (hmix(i,j) > puny) then
122 
123  ! ocean surface temperature in Kelvin
124  tsfk = sst(i,j) + tffresh
125 
126  ! shortwave radiative flux
127  swabs = (c1i - albocn) * fsw(i,j)
128 
129  ! longwave radiative flux
130  flwup = -stefan_boltzmann * tsfk**4
131 
132  ! downward latent and sensible heat fluxes
133  flh = lhcoef(i,j) * delq(i,j)
134  fsh = shcoef(i,j) * delt(i,j)
135 
136  ! compute sst change due to exchange with atm/ice above
137  ! Note: fhnet, fswthru are added in ice_therm_vertical.F
138  sst(i,j) = sst(i,j) + (c1i-aice(i,j))*dtcprho/hmix(i,j) &
139  *(fsh + flh + flwup + flw(i,j) + swabs)
140 
141  ! adjust qdp if cooling of mixed layer would occur when sst le Tf
142  if( sst(i,j) <= tf(i,j) .and. qdp(i,j) > c0i ) qdp(i,j) = c0i
143 
144  ! computed T change due to exchange with deep layers:
145  sst(i,j) = sst(i,j) - qdp(i,j)*dtcprho/hmix(i,j)
146 
147  ! compute potential to freeze or melt ice
148  frzmlt(i,j) = (tf(i,j)-sst(i,j))/dtcprho*hmix(i,j)
149  frzmlt(i,j) = min(max(frzmlt(i,j),-c1000),c1000)
150 
151  ! if sst is below freezing, reset sst to Tf
152  if (sst(i,j) <= tf(i,j)) sst(i,j) = tf(i,j)
153 
154  else
155 
156  frzmlt(i,j) = c0i
157  sst(i,j) = tf(i,j)
158 
159 ! endif ! hmix > puny
160  endif ! tmask
161  enddo ! i
162  enddo ! j
163 
164  end subroutine mixed_layer
165 
166 !=======================================================================
167 
168  end module ice_ocean
169 
170 !=======================================================================
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save lhcoef
Definition: ice_flux.f90:164
real(kind=dbl_kind), parameter c1000
subroutine mixed_layer
Definition: ice_ocean.f90:83
integer, parameter dbl_kind
real(kind=dbl_kind), dimension(:,:), allocatable, save fsw
Definition: ice_flux.f90:164
real(kind=dbl_kind), parameter c0i
logical(kind=log_kind) oceanmixed_ice
Definition: ice_ocean.f90:53
real(kind=dbl_kind), parameter stefan_boltzmann
real(kind=dbl_kind), dimension(:,:), allocatable, save flw
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save sst
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter albocn
Definition: ice_albedo.f90:48
real(kind=dbl_kind), parameter tffresh
real(kind=dbl_kind) dtice
real(kind=dbl_kind), dimension(:,:), allocatable, save hmix
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter rhow
real(kind=dbl_kind), dimension(:,:), allocatable, save qdp
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, target, save aice
Definition: ice_state.f90:82
subroutine atmo_boundary_layer(ni, sfctype, Tsf, strx, stry, Trf, Qrf, delt, delq)
Definition: ice_atmo.f90:62
real(kind=dbl_kind), dimension(:,:), allocatable, save frzmlt
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save shcoef
Definition: ice_flux.f90:164
real(kind=dbl_kind), parameter c1i
logical(kind=log_kind), dimension(:,:), allocatable tmask
Definition: ice_grid.f90:164
real(kind=dbl_kind), parameter cp_ocn