Skip to content

Commit 58f78b9

Browse files
Fill IOB, some fluxes still need to be checked
1 parent 7480492 commit 58f78b9

1 file changed

Lines changed: 53 additions & 42 deletions

File tree

config_src/mct_driver/coupler_indices.F90

Lines changed: 53 additions & 42 deletions
Original file line numberDiff line numberDiff line change
@@ -1,19 +1,19 @@
11
module coupler_indices
22

33
! MCT types
4-
use mct_mod, only : mct_aVect
4+
use mct_mod, only : mct_aVect
55
! MCT fucntions
6-
use mct_mod, only : mct_avect_indexra, mct_aVect_init, mct_aVect_clean
7-
use seq_flds_mod, only : seq_flds_x2o_fields, seq_flds_o2x_fields
8-
use seq_flds_mod, only : seq_flds_i2o_per_cat, ice_ncat
6+
use mct_mod, only : mct_avect_indexra, mct_aVect_init, mct_aVect_clean
7+
use seq_flds_mod, only : seq_flds_x2o_fields, seq_flds_o2x_fields
8+
use seq_flds_mod, only : seq_flds_i2o_per_cat, ice_ncat
99

1010
! MOM types
11-
use MOM_grid, only : ocean_grid_type
11+
use MOM_grid, only : ocean_grid_type
1212
use MOM_surface_forcing, only: ice_ocean_boundary_type
1313
! MOM functions
14-
use MOM_domains, only : pass_var, AGRID
15-
use ocean_model_mod, only : ocean_public_type
16-
14+
use MOM_domains, only : pass_var, AGRID
15+
use ocean_model_mod, only : ocean_public_type
16+
use MOM_error_handler, only : MOM_error, FATAL
1717
implicit none
1818

1919
private
@@ -116,6 +116,7 @@ subroutine coupler_indices_init(ind)
116116
call mct_aVect_init(x2o, rList=seq_flds_x2o_fields, lsize=1)
117117
call mct_aVect_init(o2x, rList=seq_flds_o2x_fields, lsize=1)
118118

119+
! ocean to coupler
119120
ind%o2x_So_t = mct_avect_indexra(o2x,'So_t')
120121
ind%o2x_So_u = mct_avect_indexra(o2x,'So_u')
121122
ind%o2x_So_v = mct_avect_indexra(o2x,'So_v')
@@ -127,14 +128,15 @@ subroutine coupler_indices_init(ind)
127128
ind%o2x_Fioo_q = mct_avect_indexra(o2x,'Fioo_q')
128129
ind%o2x_Faoo_fco2_ocn = mct_avect_indexra(o2x,'Faoo_fco2_ocn',perrWith='quiet')
129130
ind%o2x_Faoo_fdms_ocn = mct_avect_indexra(o2x,'Faoo_fdms_ocn',perrWith='quiet')
131+
132+
! coupler to ocean
130133
ind%x2o_Si_ifrac = mct_avect_indexra(x2o,'Si_ifrac')
131134
ind%x2o_Sa_pslv = mct_avect_indexra(x2o,'Sa_pslv')
132135
ind%x2o_So_duu10n = mct_avect_indexra(x2o,'So_duu10n')
133136
! QL, 150526, from wav
134137
ind%x2o_Sw_lamult = mct_avect_indexra(x2o,'Sw_lamult')
135138
ind%x2o_Sw_ustokes = mct_avect_indexra(x2o,'Sw_ustokes')
136139
ind%x2o_Sw_vstokes = mct_avect_indexra(x2o,'Sw_vstokes')
137-
138140
ind%x2o_Foxx_tauy = mct_avect_indexra(x2o,'Foxx_tauy')
139141
ind%x2o_Foxx_taux = mct_avect_indexra(x2o,'Foxx_taux')
140142
ind%x2o_Foxx_swnet = mct_avect_indexra(x2o,'Foxx_swnet')
@@ -226,7 +228,8 @@ subroutine ocn_export(ind, ocn_public, grid, o2x)
226228
do i=grid%isc,grid%iec
227229
n = n+1
228230
ig = i + grid%idg_offset
229-
o2x(ind%o2x_So_t, n) = ocn_public%t_surf(ig,jg) * grid%mask2dT(i,j)
231+
! surface temperature in Kelvin
232+
o2x(ind%o2x_So_t, n) = ocn_public%t_surf(ig,jg) * grid%mask2dT(i,j)
230233
o2x(ind%o2x_So_s, n) = ocn_public%s_surf(ig,jg) * grid%mask2dT(i,j)
231234
o2x(ind%o2x_So_u, n) = ocn_public%u_surf(ig,jg) * grid%mask2dT(i,j)
232235
o2x(ind%o2x_So_v, n) = ocn_public%v_surf(ig,jg) * grid%mask2dT(i,j)
@@ -292,15 +295,18 @@ subroutine ocn_export(ind, ocn_public, grid, o2x)
292295
end subroutine ocn_export
293296

294297

295-
subroutine fill_ice_ocean_bnd(ice_ocean_boundary, grid, x2o_o, ind)
296-
type(ice_ocean_boundary_type), intent(inout) :: ice_ocean_boundary !< A type for the ice ocean boundary
297-
type(ocean_grid_type), intent(in) :: grid
298-
!type(mct_aVect), intent(in) :: x2o_o
299-
real(kind=8), intent(in) :: x2o_o(:,:)
300-
type(cpl_indices), intent(inout) :: ind
298+
subroutine fill_ice_ocean_bnd(ice_ocean_boundary, grid, x2o_o, ind, sw_decomp, c1, c2, c3, c4)
299+
type(ice_ocean_boundary_type), intent(inout) :: ice_ocean_boundary !< A type for the ice ocean boundary
300+
type(ocean_grid_type), intent(in) :: grid
301+
!type(mct_aVect), intent(in) :: x2o_o
302+
real(kind=8), intent(in) :: x2o_o(:,:)
303+
type(cpl_indices), intent(inout) :: ind
304+
logical, intent(in) :: sw_decomp !< controls if shortwave is decomposed
305+
!! into four components
306+
real(kind=8), intent(in), optional :: c1, c2, c3, c4 !< Coeffs. used in the shortwave decomposition
301307

302308
! local variables
303-
integer :: i, j, k, ig, jg
309+
integer :: i, j, k, ig, jg !< grid indices
304310

305311
! variable that are not in ice_ocean_boundary:
306312
! latent (x2o_Foxx_lat)
@@ -324,45 +330,50 @@ subroutine fill_ice_ocean_bnd(ice_ocean_boundary, grid, x2o_o, ind)
324330

325331
! need wind_stress_multiplier?
326332

327-
! Copy from x2o to ice_ocean_boundary. ice_ocean_boundary uses global indexing with no halos.
328-
write(*,*) 'max. k is:', (grid%jec-grid%jsc+1) * (grid%iec-grid%isc+1)
329-
! zonal wind stress (taux)
330-
write(*,*) 'taux', SIZE(x2o_o(ind%x2o_Foxx_taux,:))
331-
write(*,*) 'ice_ocean_boundary%u_flux', SIZE(ice_ocean_boundary%u_flux(:,:))
332333
k = 0
333334
do j = grid%jsc, grid%jec
334335
jg = j + grid%jdg_offset
335336
do i = grid%isc, grid%iec
336337
k = k + 1 ! Increment position within gindex
337338
ig = i + grid%idg_offset
338339
! zonal wind stress (taux)
339-
ice_ocean_boundary%u_flux(ig,jg) = 0.0 ! x20_o(ind%x2o_Foxx_taux,k)
340+
ice_ocean_boundary%u_flux(ig,jg) = x2o_o(ind%x2o_Foxx_taux,k)
340341
! meridional wind stress (tauy)
341-
ice_ocean_boundary%v_flux(ig,jg) = 0.0 ! x20_o(ind%x2o_Foxx_tauy,k)
342-
! sensible heat flux
343-
ice_ocean_boundary%t_flux(ig,jg) = 0.0 ! x20_o(ind%x2o_Foxx_sen,k)
342+
ice_ocean_boundary%v_flux(ig,jg) = x2o_o(ind%x2o_Foxx_tauy,k)
343+
! sensible heat flux (W/m2)
344+
ice_ocean_boundary%t_flux(ig,jg) = -x2o_o(ind%x2o_Foxx_sen,k)
344345
! salt flux
345-
ice_ocean_boundary%salt_flux(ig,jg) = 0.0 ! x20_o(ind%x2o_Fioi_salt,k)
346+
ice_ocean_boundary%salt_flux(ig,jg) = x2o_o(ind%x2o_Fioi_salt,k)
346347
! heat flux from snow & ice melt
347-
ice_ocean_boundary%calving_hflx(ig,jg) = 0.0 ! x20_o(ind%x2o_Fioi_melth,k)
348+
ice_ocean_boundary%calving_hflx(ig,jg) = x2o_o(ind%x2o_Fioi_melth,k)
348349
! snow melt flux
349-
ice_ocean_boundary%fprec(ig,jg) = 0.0 ! x20_o(ind%x2o_Fioi_meltw,k)
350+
!ice_ocean_boundary%fprec(ig,jg) = x2o_o(ind%x2o_Fioi_meltw,k)
350351
! river runoff flux
351-
ice_ocean_boundary%runoff(ig,jg) = 0.0 ! x20_o(ind%x2o_Foxx_rofl,k)
352+
ice_ocean_boundary%runoff(ig,jg) = x2o_o(ind%x2o_Foxx_rofl,k)
352353
! ice runoff flux
353-
ice_ocean_boundary%calving(ig,jg) = 0.0 ! x20_o(ind%x2o_Foxx_rofi,k)
354+
ice_ocean_boundary%calving(ig,jg) = x2o_o(ind%x2o_Foxx_rofi,k)
354355
! liquid precipitation (rain)
355-
ice_ocean_boundary%lprec(ig,jg) = 0.0 ! x20_o(ind%x2o_Faxa_rain,k)
356-
! froze precipitation (snow)
357-
ice_ocean_boundary%fprec(ig,jg) = 0.0 ! x20_o(ind%x2o_Faxa_snow,k)
358-
!!!!!!! LONGWAVE NEEDS TO BE FIXED !!!!!!!
359-
! longwave radiation (up)
360-
ice_ocean_boundary%lw_flux(ig,jg) = 0.0 ! x20_o(k,ind%x2o_Foxx_lwup)
361-
! longwave radiation (down)
362-
ice_ocean_boundary%lw_flux(ig,jg) = 0.0 ! x20_o(k,ind%x2o_Faxa_lwdn)
363-
!!!!!!! SHORTWAVE NEEDS TO BE COMBINED !!!!!!!
364-
! net short-wave heat flux
365-
ice_ocean_boundary%u_flux(ig,jg) = 0.0 ! x20_o(k,ind%x2o_Foxx_swnet)
356+
ice_ocean_boundary%lprec(ig,jg) = x2o_o(ind%x2o_Faxa_rain,k)
357+
! frozen precipitation (snow)
358+
ice_ocean_boundary%fprec(ig,jg) = x2o_o(ind%x2o_Faxa_snow,k)
359+
! evaporation flux (kg/m2/s)
360+
ice_ocean_boundary%q_flux(ig,jg) = -x2o_o(ind%x2o_Foxx_evap,k)
361+
! longwave radiation, sum up and down (W/m2)
362+
ice_ocean_boundary%lw_flux(ig,jg) = x2o_o(ind%x2o_Faxa_lwdn,k) + x2o_o(ind%x2o_Foxx_lwup,k)
363+
if (sw_decomp) then
364+
! Use runtime coefficients to decompose net short-wave heat flux into 4 components
365+
! 1) visible, direct shortwave (W/m2)
366+
ice_ocean_boundary%sw_flux_vis_dir(ig,jg) = x2o_o(ind%x2o_Foxx_swnet,k)*c1
367+
! 2) visible, diffuse shortwave (W/m2)
368+
ice_ocean_boundary%sw_flux_vis_dif(ig,jg) = x2o_o(ind%x2o_Foxx_swnet,k)*c2
369+
! 3) near-IR, direct shortwave (W/m2)
370+
ice_ocean_boundary%sw_flux_nir_dir(ig,jg) = x2o_o(ind%x2o_Foxx_swnet,k)*c3
371+
! 4) near-IR, diffuse shortwave (W/m2)
372+
ice_ocean_boundary%sw_flux_nir_dif(ig,jg) = x2o_o(ind%x2o_Foxx_swnet,k)*c4
373+
else
374+
call MOM_error(FATAL,"fill_ice_ocean_bnd: this option has not been implemented yet."// &
375+
"Shortwave must be decomposed using coeffs. c1, c2, c3, c4.");
376+
endif
366377
enddo
367378
enddo
368379

0 commit comments

Comments
 (0)