My Project
Functions/Subroutines | Variables
ice_ocean Module Reference

Functions/Subroutines

subroutine mixed_layer
 

Variables

logical(kind=log_kind) oceanmixed_ice
 

Function/Subroutine Documentation

◆ mixed_layer()

subroutine ice_ocean::mixed_layer ( )

Definition at line 83 of file ice_ocean.f90.

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 
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
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, save fsw
Definition: ice_flux.f90:164
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), parameter stefan_boltzmann
real(kind=dbl_kind), dimension(:,:), allocatable, save flw
Definition: ice_flux.f90:91
integer(kind=int_kind) jlo
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, save sst
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter albocn
Definition: ice_albedo.f90:48
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
real(kind=dbl_kind), parameter tffresh
real(kind=dbl_kind) dtice
real(kind=dbl_kind), dimension(:,:), allocatable, save hmix
Definition: ice_flux.f90:91
integer(kind=int_kind) jhi
Definition: ice_domain.f90:101
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
Here is the call graph for this function:

Variable Documentation

◆ oceanmixed_ice

logical (kind=log_kind) ice_ocean::oceanmixed_ice

Definition at line 53 of file ice_ocean.f90.

53  logical (kind=log_kind) :: &
54  oceanmixed_ice ! if true, use ocean mixed layer