My Project
ice_flux.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_flux - flux variable declarations: coupler, diagnostic and internal
23 !
24 ! !DESCRIPTION:
25 !
26 ! Flux variable declarations; these include fields sent from the coupler
27 ! ("in"), sent to the coupler ("out"), written to diagnostic history files
28 ! ("diagnostic"), and used internally ("internal").
29 !
30 ! !REVISION HISTORY:
31 !
32 ! author Elizabeth C. Hunke, LANL
33 !
34 ! !INTERFACE:
35 !
36  module ice_flux
37 !
38 ! !USES:
39 !
40  use ice_kinds_mod
41  use ice_domain
42  use ice_constants
43 !
44 !EOP
45 !
46  implicit none
47  save
48 
49  !-----------------------------------------------------------------
50  ! Dynamics component
51  !-----------------------------------------------------------------
52 
53 ! Dynamic part using FVCOM code ! ggao
54 
55  real (kind=dbl_kind), dimension (:,:),allocatable,save :: &
56  ! real (kind=dbl_kind), dimension (ilo:ihi,jlo:jhi) :: &
57 
58  ! in from ocean
59  & ss_tltx &! sea surface slope, x-direction (m/m)
60  &, ss_tlty &! sea surface slope, y-direction
61  &, uocn &! ocean current, x-direction (m/s)
62  &, vocn &! ocean current, y-direction (m/s)
63 
64  ! out to atmos phere
65 
66  &, strairxt &! stress on ice by air, x-direction
67  &, strairyt &! stress on ice by air, y-direction
68 
69  ! out to ocean T-cell (kg/m s^2)
70  &, strocnxt &! ice-ocean stress, x-direction
71  &, strocnyt &! ice-ocean stress, y-direction
72 
73  ! diagnostic
74  &, daidtd &! &ice area tendency due to transport (s^-1)
75  &, dvidtd &! ice volume tendency due to transport (m/s)
76 
77  ! internal U-cell (kg/m s^2)
78  &, strocnx &! ice-ocean stress, x-direction
79  &, strocny &! ice-ocean stress, y-direction
80  &, strairx &! stress on ice by air, x-direction
81  &, strairy &! stress on ice by air, y-direction
82  &, strtltx &! stress due to sea surface slope, x-direction
83  &, strtlty ! stress due to sea surface slope, y-direction
84 
85  !-----------------------------------------------------------------
86  ! Thermodynamics component
87  !-----------------------------------------------------------------
88 
89  ! in from atmosphere
90 ! real (kind=dbl_kind), dimension (ilo:ihi,jlo:jhi) :: &
91  real (kind=dbl_kind), dimension (:,:),allocatable,save :: &
92 
93  & zlvl & ! atm level height (m)
94  &, uatm & ! wind speed (m/s)
95  &, vatm &
96  &, pott & ! air potential temperature (K)
97  &, tair & ! air temperature (K)
98  &, qa & ! specific humidity (kg/kg)
99  &, rhoa & ! air density (kg/m^3)
100  &, swvdr & ! sw down, visible, direct (W/m^2)
101  &, swvdf & ! sw down, visible, diffuse (W/m^2)
102  &, swidr & ! sw down, near IR, direct (W/m^2)
103  &, swidf & ! sw down, near IR, diffuse (W/m^2)
104  &, flw & ! incoming longwave radiation (W/m^2)
105  &, frain & ! rainfall rate (kg/m^2 s)
106  &, fsnow & ! snowfall rate (kg/m^2 s)
107 
108  ! in from ocean
109 
110  &, frzmlt & ! freezing/melting potential (W/m^2)
111  &, sss & ! sea surface salinity (ppt)
112  &, sst & ! sea surface temperature (C)
113  &, tf & ! freezing temperature (C)
114  &, qdp & ! deep ocean heat flux (W/m^2)
115  &, hmix ! mixed layer depth (m)
116 
117  ! out to atmosphere
118  ! note albedos are in ice_albedo.F, Tsfc in ice_state.F
119 ! real (kind=dbl_kind), dimension (ilo:ihi,jlo:jhi) :: &
120  real (kind=dbl_kind), dimension (:,:),allocatable,save :: &
121 
122  & fsens &! sensible heat flux (W/m^2)
123  &, flat &! latent heat flux (W/m^2)
124  &, fswabs &! shortwave flux absorbed in ice and ocean (W/m^2)
125  &, flwout &! outgoing longwave radiation (W/m^2)
126  &, evap &! evaporative water flux (kg/m^2/s)
127  &, tref &! 2m atm reference temperature (K)
128  &, qref &! 2m atm reference sp humidity (kg/kg)
129 
130  ! out to ocean
131 
132  &, fresh & ! fresh water flux to ocean (kg/m^2/s)
133  &, fsalt & ! salt flux to ocean (kg/m^2/s)
134  &, fhnet & ! net heat flux to ocean (W/m^2)
135  &, fswthru ! shortwave penetrating to ocean (W/m^2)
136 
137  ! diagnostic
138 
139 ! real (kind=dbl_kind), dimension (ilo:ihi,jlo:jhi) :: &
140  real (kind=dbl_kind),dimension (:,:),allocatable,save :: &
141  & congel &! basal ice growth (m/step-->cm/day)
142  &, frazil &! frazil ice growth (m/step-->cm/day)
143  &, snoice &! snow-ice formation (m/step-->cm/day)
144  &, meltt &! top ice melt (m/step-->cm/day)
145  &, meltb &! basal ice melt (m/step-->cm/day)
146  &, meltl &! lateral ice melt (m/step-->cm/day)
147  &, daidtt &! ice area tendency thermo. (s^-1)
148  &, dvidtt &! ice volume tendency thermo. (m/s)
149  &, mlt_onset &! day of year that sfc melting begins
150  &, frz_onset &! day of year that freezing begins (congel or frazil)
151 
152  ! NOTE: The following ocean diagnostic fluxes measure
153  ! the same thing as their coupler counterparts but over
154  ! different time intervals. The coupler variables are
155  ! computed from one to_coupler call to the next, whereas
156  ! the diagnostic variables are computed over a time step.
157  &, fresh_hist &! fresh water flux to ocean (kg/m^2/s)
158  &, fsalt_hist &! salt flux to ocean (kg/m^2/s)
159  &, fhnet_hist &! net heat flux to ocean (W/m^2)
160  &, fswthru_hist ! shortwave penetrating to ocean (W/m^2)
161 
162  ! internal
163 ! real (kind=dbl_kind), dimension (ilo:ihi,jlo:jhi) :: &
164  real (kind=dbl_kind),dimension (:,:),allocatable,save :: &
165  & fsw &! incoming shortwave radiation (W/m^2)
166  &, wind &! wind speed (m/s)
167  &, shcoef &! transfer coefficient for sensible heat
168  &, lhcoef ! transfer coefficient for latent heat
169 
170 !=======================================================================
171 
172  contains
173 
174 !=======================================================================
175 !BOP
176 !
177 ! !IROUTINE: init_flux_atm - initialize all atmospheric fluxes sent to coupler
178 !
179 ! !INTERFACE:
180 !
181  subroutine init_flux_atm
182 !
183 ! !DESCRIPTION:
184 !
185 ! Initialize all fluxes sent to coupler for use by the atm model
186 ! and a few state quantities
187 !
188 ! !REVISION HISTORY:
189 !
190 ! author: Elizabeth C. Hunke, LANL
191 !
192 ! !USES:
193  use ice_state, only : aice
194 !
195 ! !INPUT/OUTPUT PARAMETERS:
196 !
197 !
198 !EOP
199 !
200  integer (kind=int_kind) :: &
201  & i, j ! horizontal indices
202 
203  !-----------------------------------------------------------------
204  ! fluxes sent
205  !-----------------------------------------------------------------
206 
207  do j = jlo, jhi
208  do i = ilo, ihi
209 
210  strairxt(i,j) = c0i ! wind stress, T grid
211  strairyt(i,j) = c0i
212  fsens(i,j) = c0i
213  flat(i,j) = c0i
214  fswabs(i,j) = c0i
215  flwout(i,j) = c0i
216  evap(i,j) = c0i
217  tref(i,j) = c0i
218  qref(i,j) = c0i
219 
220  !-----------------------------------------------------------------
221  ! other miscellaneous fields
222  !-----------------------------------------------------------------
223  strairx(i,j) = c0i ! wind stress, U grid
224  strairy(i,j) = c0i
225 
226  enddo
227  enddo
228 
229  end subroutine init_flux_atm
230 
231 !=======================================================================
232 !BOP
233 !
234 ! !IROUTINE: init_flux_ocn - initialize ocean fluxes sent to coupler
235 !
236 ! !INTERFACE:
237 !
238  subroutine init_flux_ocn
240 !cdir$ inlinenever init_flux_ocn
241 !
242 ! !DESCRIPTION:
243 !
244 ! Initialize fluxes sent to coupler for use by the ocean model
245 !
246 ! !REVISION HISTORY:
247 !
248 ! author: Elizabeth C. Hunke, LANL
249 !
250 ! !USES:
251 !
252 ! !INPUT/OUTPUT PARAMETERS:
253 !
254 !
255 !EOP
256 !
257  integer (kind=int_kind) :: &
258  & i, j ! horizontal indices
259 
260  !-----------------------------------------------------------------
261  ! fluxes sent
262  !-----------------------------------------------------------------
263 
264  do j = jlo, jhi
265  do i = ilo, ihi
266  fresh(i,j) = c0i
267  fsalt(i,j) = c0i
268  fhnet(i,j) = c0i
269  fswthru(i,j) = c0i
270 
271  qdp(i,j) = c0i
272  hmix(i,j) = c20
273  enddo
274  enddo
275 
276  end subroutine init_flux_ocn
277 
278 !=======================================================================
279 !BOP
280 !
281 ! !IROUTINE: init_diagnostics - initialize diagnostic fields
282 !
283 ! !INTERFACE:
284 !
285  subroutine init_diagnostics
286 !
287 ! !DESCRIPTION:
288 !
289 ! Initialize diagnostic fields written to history files.
290 !
291 ! !REVISION HISTORY:
292 !
293 ! author: William H. Lipscomb, LANL
294 !
295 ! !USES:
296 !
297 ! !INPUT/OUTPUT PARAMETERS:
298 !
299 !
300 !EOP
301 !
302  use ice_state, only: aice, vice
303 
304  integer (kind=int_kind) :: &
305  & i, j ! horizontal indices
306 
307  do j = jlo, jhi
308  do i = ilo, ihi
309  congel(i,j) = c0i
310  frazil(i,j) = c0i
311  snoice(i,j) = c0i
312  meltt(i,j) = c0i
313  meltb(i,j) = c0i
314  meltl(i,j) = c0i
315  daidtt(i,j) = aice(i,j) ! temporarily used for initial area
316  dvidtt(i,j) = vice(i,j) ! temporarily used for initial volume
317  daidtd(i,j) = c0i
318  dvidtd(i,j) = c0i
319  fresh_hist(i,j) = c0i
320  fsalt_hist(i,j) = c0i
321  fhnet_hist(i,j) = c0i
322  fswthru_hist(i,j) = c0i
323  enddo ! i
324  enddo ! j
325 
326  end subroutine init_diagnostics
327 
328 !=======================================================================
329 !BOP
330 !
331 ! !IROUTINE: merge_fluxes - aggregate flux information over ITD
332 !
333 ! !INTERFACE:
334 !
335  subroutine merge_fluxes (ni, &
336  & strxn, stryn, fsensn, flatn, fswabsn, &
337  & flwoutn, evapn, Trefn, Qrefn, &
338  & freshn, fsaltn, fhnetn, fswthrun)
339 !
340 ! !DESCRIPTION:
341 !
342 ! Aggregates flux information from all ice thickness categories
343 !
344 ! !REVISION HISTORY:
345 !
346 ! author: Elizabeth C. Hunke, LANL
347 !
348 ! !USES:
349 !
350  use ice_state
351 !
352 ! !INPUT/OUTPUT PARAMETERS:
353 !
354  integer (kind=int_kind), intent(in) :: &
355  & ni ! thickness category index
356 
357  real (kind=dbl_kind), dimension(ilo:ihi,jlo:jhi), intent(in) :: &
358  & strxn &! air/ice zonal strss, (N/m**2)
359  &, stryn &! air/ice merdnl strss, (N/m**2)
360  &, fsensn &! sensible heat flx (W/m**2)
361  &, flatn &! latent heat flx (W/m**2)
362  &, fswabsn &! shortwave absorbed heat flx (W/m**2)
363  &, flwoutn &! upwd lw emitted heat flx (W/m**2)
364  &, evapn &! evaporation (kg/m2/s)
365  &, trefn &! air tmp reference level (K)
366  &, qrefn &! air sp hum reference level (kg/kg)
367  &, freshn &! fresh water flux to ocean (kg/m2/s)
368  &, fsaltn &! salt flux to ocean (kg/m2/s)
369  &, fhnetn &! actual ocn/ice heat flx (W/m**2)
370  &, fswthrun ! sw radiation through ice bot (W/m**2)
371 !
372 !EOP
373 !
374  integer (kind=int_kind) :: &
375  & i, j ! horizontal indices
376 
377  do j = jlo,jhi
378  do i = ilo,ihi
379 
380  ! atmo fluxes
381 
382  strairxt(i,j) = strairxt(i,j) + strxn(i,j) * aicen(i,j,ni)
383  strairyt(i,j) = strairyt(i,j) + stryn(i,j) * aicen(i,j,ni)
384  fsens(i,j) = fsens(i,j) + fsensn(i,j) * aicen(i,j,ni)
385  flat(i,j) = flat(i,j) + flatn(i,j) * aicen(i,j,ni)
386  fswabs(i,j) = fswabs(i,j) + fswabsn(i,j) * aicen(i,j,ni)
387  flwout(i,j) = flwout(i,j) + &
388  & (flwoutn(i,j) - (c1i-emissivity)*flw(i,j))* aicen(i,j,ni)
389  evap(i,j) = evap(i,j) + evapn(i,j) * aicen(i,j,ni)
390  tref(i,j) = tref(i,j) + trefn(i,j) * aicen(i,j,ni)
391  qref(i,j) = qref(i,j) + qrefn(i,j) * aicen(i,j,ni)
392 
393  ! ocean fluxes: update both coupler and history variables
394 
395  fresh(i,j) = fresh(i,j) + freshn(i,j) * aicen(i,j,ni)
396  fresh_hist(i,j) = fresh_hist(i,j) + freshn(i,j) * aicen(i,j,ni)
397  fsalt(i,j) = fsalt(i,j) + fsaltn(i,j) * aicen(i,j,ni)
398  fsalt_hist(i,j) = fsalt_hist(i,j) + fsaltn(i,j) * aicen(i,j,ni)
399  fhnet(i,j) = fhnet(i,j) + fhnetn(i,j) * aicen(i,j,ni)
400  fhnet_hist(i,j) = fhnet_hist(i,j) + fhnetn(i,j) * aicen(i,j,ni)
401  fswthru(i,j) = fswthru(i,j) + fswthrun(i,j) * aicen(i,j,ni)
402  fswthru_hist(i,j) = fswthru_hist(i,j) &
403  & + fswthrun(i,j) * aicen(i,j,ni)
404  enddo ! i
405  enddo ! j
406 
407  end subroutine merge_fluxes
408 
409 !=======================================================================
410 
411  end module ice_flux
412 
413 !=======================================================================
real(kind=dbl_kind), dimension(:,:), allocatable, save strairyt
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save tf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save fsalt_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save lhcoef
Definition: ice_flux.f90:164
real(kind=dbl_kind), dimension(:,:), allocatable, save strtltx
Definition: ice_flux.f90:55
integer, parameter dbl_kind
real(kind=dbl_kind), dimension(:,:), allocatable, save snoice
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save strocnx
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save strocnxt
Definition: ice_flux.f90:55
integer(kind=int_kind) ihi
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, save sss
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save fsw
Definition: ice_flux.f90:164
real(kind=dbl_kind), parameter c0i
real(kind=dbl_kind), dimension(:,:), allocatable, save uocn
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save meltl
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save swvdf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save fhnet_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save vocn
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save rhoa
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save pott
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save flwout
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save fsens
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save fresh
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save flw
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter c20
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), dimension(:,:), allocatable, save congel
Definition: ice_flux.f90:140
subroutine init_flux_atm
Definition: ice_flux.f90:182
real(kind=dbl_kind), dimension(:,:), allocatable, save vatm
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter emissivity
real(kind=dbl_kind), dimension(:,:), allocatable, save daidtd
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save ss_tltx
Definition: ice_flux.f90:55
integer(kind=int_kind) ilo
Definition: ice_domain.f90:101
real(kind=dbl_kind), dimension(:,:), allocatable, save ss_tlty
Definition: ice_flux.f90:55
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), dimension(:,:), allocatable, save meltt
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save qref
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save fswthru
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save frain
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save frazil
Definition: ice_flux.f90:140
subroutine merge_fluxes(ni, strxn, stryn, fsensn, flatn, fswabsn, flwoutn, evapn, Trefn, Qrefn, freshn, fsaltn, fhnetn, fswthrun)
Definition: ice_flux.f90:339
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
real(kind=dbl_kind), dimension(:,:), allocatable, save fresh_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save fsalt
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save flat
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save fswabs
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save strairy
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save strairxt
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:,:), allocatable, target, save aicen
Definition: ice_state.f90:97
real(kind=dbl_kind), dimension(:,:), allocatable, save fswthru_hist
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save dvidtt
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save evap
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save uatm
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save daidtt
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save swvdr
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save strtlty
Definition: ice_flux.f90:55
subroutine init_diagnostics
Definition: ice_flux.f90:286
real(kind=dbl_kind), dimension(:,:), allocatable, save frzmlt
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save zlvl
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save meltb
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save shcoef
Definition: ice_flux.f90:164
real(kind=dbl_kind), dimension(:,:), allocatable, target, save vice
Definition: ice_state.f90:82
real(kind=dbl_kind), dimension(:,:), allocatable, save strocnyt
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save fsnow
Definition: ice_flux.f90:91
real(kind=dbl_kind), parameter c1i
real(kind=dbl_kind), dimension(:,:), allocatable, save strocny
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save qa
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save swidf
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save frz_onset
Definition: ice_flux.f90:140
subroutine init_flux_ocn
Definition: ice_flux.f90:239
real(kind=dbl_kind), dimension(:,:), allocatable, save tref
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save wind
Definition: ice_flux.f90:164
real(kind=dbl_kind), dimension(:,:), allocatable, save tair
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save swidr
Definition: ice_flux.f90:91
real(kind=dbl_kind), dimension(:,:), allocatable, save mlt_onset
Definition: ice_flux.f90:140
real(kind=dbl_kind), dimension(:,:), allocatable, save dvidtd
Definition: ice_flux.f90:55
real(kind=dbl_kind), dimension(:,:), allocatable, save fhnet
Definition: ice_flux.f90:120
real(kind=dbl_kind), dimension(:,:), allocatable, save strairx
Definition: ice_flux.f90:55