1 module mom_surface_forcing_gfdl
10 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
14 use mom_domains,
only : agrid, bgrid_ne, cgrid_ne, to_all
15 use mom_domains,
only : to_north, to_east, omit_corners
26 use mom_restart,
only : restart_init_end, save_restart, restore_state
34 use coupler_types_mod,
only : coupler_2d_bc_type, coupler_type_write_chksums
35 use coupler_types_mod,
only : coupler_type_initialized, coupler_type_spawn
36 use coupler_types_mod,
only : coupler_type_copy_data
37 use data_override_mod,
only : data_override_init, data_override
38 use fms_mod,
only : stdout
39 use mpp_mod,
only : mpp_chksum
40 use time_interp_external_mod,
only : init_external_field, time_interp_external
41 use time_interp_external_mod,
only : time_interp_external_init
43 implicit none ;
private
45 #include <MOM_memory.h>
47 public convert_iob_to_fluxes, convert_iob_to_forces
48 public surface_forcing_init
49 public ice_ocn_bnd_type_chksum
50 public forcing_save_restart
60 integer :: wind_stagger
64 logical :: use_temperature
65 real :: wind_stress_multiplier
68 real :: area_surf = -1.0
69 real :: latent_heat_fusion
70 real :: latent_heat_vapor
78 logical :: use_limited_p_ssh
82 logical :: approx_net_mass_src
87 logical :: read_gust_2d
88 real,
pointer,
dimension(:,:) :: &
91 real,
pointer,
dimension(:,:) :: &
94 real,
pointer,
dimension(:,:) :: &
98 logical :: read_tideamp
100 logical :: rigid_sea_ice
104 real :: density_sea_ice
106 real :: rigid_sea_ice_mass
108 logical :: allow_flux_adjustments
110 logical :: restore_salt
112 logical :: restore_temp
115 logical :: salt_restore_as_sflux
116 logical :: adjust_net_srestore_to_zero
117 logical :: adjust_net_srestore_by_scaling
118 logical :: adjust_net_fresh_water_to_zero
119 logical :: use_net_fw_adjustment_sign_bug
120 logical :: adjust_net_fresh_water_by_scaling
121 logical :: mask_srestore_under_ice
123 real :: ice_salt_concentration
124 logical :: mask_srestore_marginal_seas
125 real :: max_delta_srestore
126 real :: max_delta_trestore
127 real,
pointer,
dimension(:,:) :: basin_mask => null()
128 logical :: answers_2018
131 logical :: fix_ustar_gustless_bug
133 logical :: check_no_land_fluxes
136 character(len=200) :: inputdir
137 character(len=200) :: salt_restore_file
138 character(len=30) :: salt_restore_var_name
139 logical :: mask_srestore
143 real,
pointer,
dimension(:,:) :: srestore_mask => null()
144 character(len=200) :: temp_restore_file
145 character(len=30) :: temp_restore_var_name
146 logical :: mask_trestore
150 real,
pointer,
dimension(:,:) :: trestore_mask => null()
151 integer :: id_srestore = -1
152 integer :: id_trestore = -1
165 real,
pointer,
dimension(:,:) :: u_flux =>null()
166 real,
pointer,
dimension(:,:) :: v_flux =>null()
167 real,
pointer,
dimension(:,:) :: t_flux =>null()
168 real,
pointer,
dimension(:,:) :: q_flux =>null()
169 real,
pointer,
dimension(:,:) :: salt_flux =>null()
170 real,
pointer,
dimension(:,:) :: lw_flux =>null()
171 real,
pointer,
dimension(:,:) :: sw_flux_vis_dir =>null()
172 real,
pointer,
dimension(:,:) :: sw_flux_vis_dif =>null()
173 real,
pointer,
dimension(:,:) :: sw_flux_nir_dir =>null()
174 real,
pointer,
dimension(:,:) :: sw_flux_nir_dif =>null()
175 real,
pointer,
dimension(:,:) :: lprec =>null()
176 real,
pointer,
dimension(:,:) :: fprec =>null()
177 real,
pointer,
dimension(:,:) :: runoff =>null()
178 real,
pointer,
dimension(:,:) :: calving =>null()
179 real,
pointer,
dimension(:,:) :: stress_mag =>null()
180 real,
pointer,
dimension(:,:) :: ustar_berg =>null()
181 real,
pointer,
dimension(:,:) :: area_berg =>null()
182 real,
pointer,
dimension(:,:) :: mass_berg =>null()
183 real,
pointer,
dimension(:,:) :: runoff_hflx =>null()
184 real,
pointer,
dimension(:,:) :: calving_hflx =>null()
185 real,
pointer,
dimension(:,:) :: p =>null()
187 real,
pointer,
dimension(:,:) :: mi =>null()
188 real,
pointer,
dimension(:,:) :: ice_rigidity =>null()
193 type(coupler_2d_bc_type) :: fluxes
195 integer :: wind_stagger = -999
201 integer :: id_clock_forcing
208 subroutine convert_iob_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G, US, CS, sfc_state)
210 target,
intent(in) :: iob
212 type(
forcing),
intent(inout) :: fluxes
215 integer,
dimension(4),
intent(in) :: index_bounds
216 type(time_type),
intent(in) :: time
218 real,
intent(in) :: valid_time
224 type(
surface),
intent(in) :: sfc_state
227 real,
dimension(SZI_(G),SZJ_(G)) :: &
240 integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq, i0, j0
241 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, isr, ier, jsr, jer
242 integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
247 real :: kg_m2_s_conversion
251 real :: sign_for_net_fw_bug
253 call cpu_clock_begin(id_clock_forcing)
255 isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
256 jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
257 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
258 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
259 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
260 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
261 isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
263 kg_m2_s_conversion = us%kg_m2s_to_RZ_T
264 if (cs%restore_temp) rhoxcp = cs%Rho0 * fluxes%C_p
265 open_ocn_mask(:,:) = 1.0
267 fluxes%vPrecGlobalAdj = 0.0
268 fluxes%vPrecGlobalScl = 0.0
269 fluxes%saltFluxGlobalAdj = 0.0
270 fluxes%saltFluxGlobalScl = 0.0
271 fluxes%netFWGlobalAdj = 0.0
272 fluxes%netFWGlobalScl = 0.0
276 if (fluxes%dt_buoy_accum < 0)
then
278 fix_accum_bug=cs%fix_ustar_gustless_bug)
280 call safe_alloc_ptr(fluxes%sw_vis_dir,isd,ied,jsd,jed)
281 call safe_alloc_ptr(fluxes%sw_vis_dif,isd,ied,jsd,jed)
282 call safe_alloc_ptr(fluxes%sw_nir_dir,isd,ied,jsd,jed)
283 call safe_alloc_ptr(fluxes%sw_nir_dif,isd,ied,jsd,jed)
285 call safe_alloc_ptr(fluxes%p_surf,isd,ied,jsd,jed)
286 call safe_alloc_ptr(fluxes%p_surf_full,isd,ied,jsd,jed)
287 if (cs%use_limited_P_SSH)
then
288 fluxes%p_surf_SSH => fluxes%p_surf
290 fluxes%p_surf_SSH => fluxes%p_surf_full
293 call safe_alloc_ptr(fluxes%salt_flux,isd,ied,jsd,jed)
294 call safe_alloc_ptr(fluxes%salt_flux_in,isd,ied,jsd,jed)
295 call safe_alloc_ptr(fluxes%salt_flux_added,isd,ied,jsd,jed)
297 call safe_alloc_ptr(fluxes%TKE_tidal,isd,ied,jsd,jed)
298 call safe_alloc_ptr(fluxes%ustar_tidal,isd,ied,jsd,jed)
300 if (cs%allow_flux_adjustments)
then
301 call safe_alloc_ptr(fluxes%heat_added,isd,ied,jsd,jed)
302 call safe_alloc_ptr(fluxes%salt_flux_added,isd,ied,jsd,jed)
305 do j=js-2,je+2 ;
do i=is-2,ie+2
306 fluxes%TKE_tidal(i,j) = cs%TKE_tidal(i,j)
307 fluxes%ustar_tidal(i,j) = cs%ustar_tidal(i,j)
310 if (cs%restore_temp)
call safe_alloc_ptr(fluxes%heat_added,isd,ied,jsd,jed)
315 if (((
associated(iob%ustar_berg) .and. (.not.
associated(fluxes%ustar_berg))) &
316 .or. (
associated(iob%area_berg) .and. (.not.
associated(fluxes%area_berg)))) &
317 .or. (
associated(iob%mass_berg) .and. (.not.
associated(fluxes%mass_berg)))) &
320 if ((.not.coupler_type_initialized(fluxes%tr_fluxes)) .and. &
321 coupler_type_initialized(iob%fluxes)) &
322 call coupler_type_spawn(iob%fluxes, fluxes%tr_fluxes, &
323 (/is,is,ie,ie/), (/js,js,je,je/))
329 if (cs%area_surf < 0.0)
then
330 do j=js,je ;
do i=is,ie
331 work_sum(i,j) = us%L_to_m**2*g%areaT(i,j) * g%mask2dT(i,j)
338 fluxes%fluxes_used = .false.
339 fluxes%dt_buoy_accum = us%s_to_T*valid_time
341 if (cs%allow_flux_adjustments)
then
342 fluxes%heat_added(:,:) = 0.0
343 fluxes%salt_flux_added(:,:) = 0.0
346 do j=js,je ;
do i=is,ie
347 fluxes%salt_flux(i,j) = 0.0
348 fluxes%vprec(i,j) = 0.0
352 if (cs%restore_salt)
then
353 call time_interp_external(cs%id_srestore,time,data_restore)
355 open_ocn_mask(:,:) = 1.0
356 if (cs%mask_srestore_under_ice)
then
357 do j=js,je ;
do i=is,ie
358 if (sfc_state%SST(i,j) <= -0.0539*sfc_state%SSS(i,j)) open_ocn_mask(i,j)=0.0
361 if (cs%salt_restore_as_sflux)
then
362 do j=js,je ;
do i=is,ie
363 delta_sss = data_restore(i,j)- sfc_state%SSS(i,j)
364 delta_sss = sign(1.0,delta_sss)*min(abs(delta_sss),cs%max_delta_srestore)
365 fluxes%salt_flux(i,j) = 1.e-3*g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)* &
366 (cs%basin_mask(i,j)*open_ocn_mask(i,j)*cs%srestore_mask(i,j)) *delta_sss
368 if (cs%adjust_net_srestore_to_zero)
then
369 if (cs%adjust_net_srestore_by_scaling)
then
370 call adjust_area_mean_to_zero(fluxes%salt_flux, g, fluxes%saltFluxGlobalScl, &
371 unit_scale=us%RZ_T_to_kg_m2s)
372 fluxes%saltFluxGlobalAdj = 0.
374 work_sum(is:ie,js:je) = us%L_to_m**2*us%RZ_T_to_kg_m2s * &
375 g%areaT(is:ie,js:je)*fluxes%salt_flux(is:ie,js:je)
376 fluxes%saltFluxGlobalAdj =
reproducing_sum(work_sum(:,:), isr,ier, jsr,jer)/cs%area_surf
377 fluxes%salt_flux(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je) - kg_m2_s_conversion * fluxes%saltFluxGlobalAdj
380 fluxes%salt_flux_added(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je)
382 do j=js,je ;
do i=is,ie
383 if (g%mask2dT(i,j) > 0.5)
then
384 delta_sss = sfc_state%SSS(i,j) - data_restore(i,j)
385 delta_sss = sign(1.0,delta_sss)*min(abs(delta_sss),cs%max_delta_srestore)
386 fluxes%vprec(i,j) = (cs%basin_mask(i,j)*open_ocn_mask(i,j)*cs%srestore_mask(i,j))* &
387 (cs%Rho0*cs%Flux_const) * &
388 delta_sss / (0.5*(sfc_state%SSS(i,j) + data_restore(i,j)))
391 if (cs%adjust_net_srestore_to_zero)
then
392 if (cs%adjust_net_srestore_by_scaling)
then
393 call adjust_area_mean_to_zero(fluxes%vprec, g, fluxes%vPrecGlobalScl, &
394 unit_scale=us%RZ_T_to_kg_m2s)
395 fluxes%vPrecGlobalAdj = 0.
397 work_sum(is:ie,js:je) = us%L_to_m**2*g%areaT(is:ie,js:je) * &
398 us%RZ_T_to_kg_m2s*fluxes%vprec(is:ie,js:je)
399 fluxes%vPrecGlobalAdj =
reproducing_sum(work_sum(:,:), isr, ier, jsr, jer) / cs%area_surf
400 do j=js,je ;
do i=is,ie
401 fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion*fluxes%vPrecGlobalAdj ) * g%mask2dT(i,j)
409 if (cs%restore_temp)
then
410 call time_interp_external(cs%id_trestore,time,data_restore)
411 do j=js,je ;
do i=is,ie
412 delta_sst = data_restore(i,j)- sfc_state%SST(i,j)
413 delta_sst = sign(1.0,delta_sst)*min(abs(delta_sst),cs%max_delta_trestore)
414 fluxes%heat_added(i,j) = g%mask2dT(i,j) * cs%trestore_mask(i,j) * &
415 rhoxcp * delta_sst * cs%Flux_const
421 i0 = is - isc_bnd ; j0 = js - jsc_bnd
422 do j=js,je ;
do i=is,ie
424 if (
associated(iob%lprec))
then
425 fluxes%lprec(i,j) = kg_m2_s_conversion * iob%lprec(i-i0,j-j0) * g%mask2dT(i,j)
426 if (cs%check_no_land_fluxes) &
427 call check_mask_val_consistency(iob%lprec(i-i0,j-j0), g%mask2dT(i,j), i, j,
'lprec', g)
430 if (
associated(iob%fprec))
then
431 fluxes%fprec(i,j) = kg_m2_s_conversion * iob%fprec(i-i0,j-j0) * g%mask2dT(i,j)
432 if (cs%check_no_land_fluxes) &
433 call check_mask_val_consistency(iob%fprec(i-i0,j-j0), g%mask2dT(i,j), i, j,
'fprec', g)
436 if (
associated(iob%q_flux))
then
437 fluxes%evap(i,j) = - kg_m2_s_conversion * iob%q_flux(i-i0,j-j0) * g%mask2dT(i,j)
438 if (cs%check_no_land_fluxes) &
439 call check_mask_val_consistency(iob%q_flux(i-i0,j-j0), g%mask2dT(i,j), i, j,
'q_flux', g)
442 if (
associated(iob%runoff))
then
443 fluxes%lrunoff(i,j) = kg_m2_s_conversion * iob%runoff(i-i0,j-j0) * g%mask2dT(i,j)
444 if (cs%check_no_land_fluxes) &
445 call check_mask_val_consistency(iob%runoff(i-i0,j-j0), g%mask2dT(i,j), i, j,
'runoff', g)
448 if (
associated(iob%calving))
then
449 fluxes%frunoff(i,j) = kg_m2_s_conversion * iob%calving(i-i0,j-j0) * g%mask2dT(i,j)
450 if (cs%check_no_land_fluxes) &
451 call check_mask_val_consistency(iob%calving(i-i0,j-j0), g%mask2dT(i,j), i, j,
'calving', g)
454 if (
associated(iob%ustar_berg))
then
455 fluxes%ustar_berg(i,j) = us%m_to_Z*us%T_to_s * iob%ustar_berg(i-i0,j-j0) * g%mask2dT(i,j)
456 if (cs%check_no_land_fluxes) &
457 call check_mask_val_consistency(iob%ustar_berg(i-i0,j-j0), g%mask2dT(i,j), i, j,
'ustar_berg', g)
460 if (
associated(iob%area_berg))
then
461 fluxes%area_berg(i,j) = iob%area_berg(i-i0,j-j0) * g%mask2dT(i,j)
462 if (cs%check_no_land_fluxes) &
463 call check_mask_val_consistency(iob%area_berg(i-i0,j-j0), g%mask2dT(i,j), i, j,
'area_berg', g)
466 if (
associated(iob%mass_berg))
then
467 fluxes%mass_berg(i,j) = us%m_to_Z*us%kg_m3_to_R * iob%mass_berg(i-i0,j-j0) * g%mask2dT(i,j)
468 if (cs%check_no_land_fluxes) &
469 call check_mask_val_consistency(iob%mass_berg(i-i0,j-j0), g%mask2dT(i,j), i, j,
'mass_berg', g)
472 if (
associated(iob%runoff_hflx))
then
473 fluxes%heat_content_lrunoff(i,j) = us%W_m2_to_QRZ_T * iob%runoff_hflx(i-i0,j-j0) * g%mask2dT(i,j)
474 if (cs%check_no_land_fluxes) &
475 call check_mask_val_consistency(iob%runoff_hflx(i-i0,j-j0), g%mask2dT(i,j), i, j,
'runoff_hflx', g)
478 if (
associated(iob%calving_hflx))
then
479 fluxes%heat_content_frunoff(i,j) = us%W_m2_to_QRZ_T * iob%calving_hflx(i-i0,j-j0) * g%mask2dT(i,j)
480 if (cs%check_no_land_fluxes) &
481 call check_mask_val_consistency(iob%calving_hflx(i-i0,j-j0), g%mask2dT(i,j), i, j,
'calving_hflx', g)
484 if (
associated(iob%lw_flux))
then
485 fluxes%LW(i,j) = us%W_m2_to_QRZ_T * iob%lw_flux(i-i0,j-j0) * g%mask2dT(i,j)
486 if (cs%check_no_land_fluxes) &
487 call check_mask_val_consistency(iob%lw_flux(i-i0,j-j0), g%mask2dT(i,j), i, j,
'lw_flux', g)
490 if (
associated(iob%t_flux))
then
491 fluxes%sens(i,j) = -us%W_m2_to_QRZ_T* iob%t_flux(i-i0,j-j0) * g%mask2dT(i,j)
492 if (cs%check_no_land_fluxes) &
493 call check_mask_val_consistency(iob%t_flux(i-i0,j-j0), g%mask2dT(i,j), i, j,
't_flux', g)
496 fluxes%latent(i,j) = 0.0
497 if (
associated(iob%fprec))
then
498 fluxes%latent(i,j) = fluxes%latent(i,j) - &
499 iob%fprec(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_fusion
500 fluxes%latent_fprec_diag(i,j) = -g%mask2dT(i,j) * iob%fprec(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_fusion
502 if (
associated(iob%calving))
then
503 fluxes%latent(i,j) = fluxes%latent(i,j) - &
504 iob%calving(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_fusion
505 fluxes%latent_frunoff_diag(i,j) = -g%mask2dT(i,j) * iob%calving(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_fusion
507 if (
associated(iob%q_flux))
then
508 fluxes%latent(i,j) = fluxes%latent(i,j) - &
509 iob%q_flux(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_vapor
510 fluxes%latent_evap_diag(i,j) = -g%mask2dT(i,j) * iob%q_flux(i-i0,j-j0)*us%W_m2_to_QRZ_T*cs%latent_heat_vapor
513 fluxes%latent(i,j) = g%mask2dT(i,j) * fluxes%latent(i,j)
515 if (
associated(iob%sw_flux_vis_dir))
then
516 fluxes%sw_vis_dir(i,j) = g%mask2dT(i,j) * us%W_m2_to_QRZ_T * iob%sw_flux_vis_dir(i-i0,j-j0)
517 if (cs%check_no_land_fluxes) &
518 call check_mask_val_consistency(iob%sw_flux_vis_dir(i-i0,j-j0), g%mask2dT(i,j), i, j,
'sw_flux_vis_dir', g)
520 if (
associated(iob%sw_flux_vis_dif))
then
521 fluxes%sw_vis_dif(i,j) = g%mask2dT(i,j) * us%W_m2_to_QRZ_T * iob%sw_flux_vis_dif(i-i0,j-j0)
522 if (cs%check_no_land_fluxes) &
523 call check_mask_val_consistency(iob%sw_flux_vis_dif(i-i0,j-j0), g%mask2dT(i,j), i, j,
'sw_flux_vis_dif', g)
525 if (
associated(iob%sw_flux_nir_dir))
then
526 fluxes%sw_nir_dir(i,j) = g%mask2dT(i,j) * us%W_m2_to_QRZ_T * iob%sw_flux_nir_dir(i-i0,j-j0)
527 if (cs%check_no_land_fluxes) &
528 call check_mask_val_consistency(iob%sw_flux_nir_dir(i-i0,j-j0), g%mask2dT(i,j), i, j,
'sw_flux_nir_dir', g)
530 if (
associated(iob%sw_flux_nir_dif))
then
531 fluxes%sw_nir_dif(i,j) = g%mask2dT(i,j) * us%W_m2_to_QRZ_T * iob%sw_flux_nir_dif(i-i0,j-j0)
532 if (cs%check_no_land_fluxes) &
533 call check_mask_val_consistency(iob%sw_flux_nir_dif(i-i0,j-j0), g%mask2dT(i,j), i, j,
'sw_flux_nir_dif', g)
535 if (cs%answers_2018)
then
536 fluxes%sw(i,j) = fluxes%sw_vis_dir(i,j) + fluxes%sw_vis_dif(i,j) + &
537 fluxes%sw_nir_dir(i,j) + fluxes%sw_nir_dif(i,j)
539 fluxes%sw(i,j) = (fluxes%sw_vis_dir(i,j) + fluxes%sw_vis_dif(i,j)) + &
540 (fluxes%sw_nir_dir(i,j) + fluxes%sw_nir_dif(i,j))
546 if (
associated(iob%p))
then
547 if (cs%max_p_surf >= 0.0)
then
548 do j=js,je ;
do i=is,ie
549 fluxes%p_surf_full(i,j) = g%mask2dT(i,j) * us%kg_m3_to_R*us%m_s_to_L_T**2*iob%p(i-i0,j-j0)
550 fluxes%p_surf(i,j) = min(fluxes%p_surf_full(i,j),cs%max_p_surf)
551 if (cs%check_no_land_fluxes) &
552 call check_mask_val_consistency(iob%p(i-i0,j-j0), g%mask2dT(i,j), i, j,
'p', g)
555 do j=js,je ;
do i=is,ie
556 fluxes%p_surf_full(i,j) = g%mask2dT(i,j) * us%kg_m3_to_R*us%m_s_to_L_T**2*iob%p(i-i0,j-j0)
557 fluxes%p_surf(i,j) = fluxes%p_surf_full(i,j)
558 if (cs%check_no_land_fluxes) &
559 call check_mask_val_consistency(iob%p(i-i0,j-j0), g%mask2dT(i,j), i, j,
'p', g)
562 fluxes%accumulate_p_surf = .true.
566 if (
associated(iob%salt_flux))
then
567 do j=js,je ;
do i=is,ie
568 fluxes%salt_flux(i,j) = g%mask2dT(i,j)*(fluxes%salt_flux(i,j) - kg_m2_s_conversion*iob%salt_flux(i-i0,j-j0))
569 fluxes%salt_flux_in(i,j) = g%mask2dT(i,j)*( -kg_m2_s_conversion*iob%salt_flux(i-i0,j-j0) )
570 if (cs%check_no_land_fluxes) &
571 call check_mask_val_consistency(iob%salt_flux(i-i0,j-j0), g%mask2dT(i,j), i, j,
'salt_flux', g)
586 if (cs%adjust_net_fresh_water_to_zero)
then
587 sign_for_net_fw_bug = 1.
588 if (cs%use_net_FW_adjustment_sign_bug) sign_for_net_fw_bug = -1.
589 do j=js,je ;
do i=is,ie
590 net_fw(i,j) = us%RZ_T_to_kg_m2s* &
591 (((fluxes%lprec(i,j) + fluxes%fprec(i,j)) + &
592 (fluxes%lrunoff(i,j) + fluxes%frunoff(i,j))) + &
593 (fluxes%evap(i,j) + fluxes%vprec(i,j)) ) * us%L_to_m**2*g%areaT(i,j)
600 if (
associated(iob%salt_flux) .and. (cs%ice_salt_concentration>0.0)) &
601 net_fw(i,j) = net_fw(i,j) + sign_for_net_fw_bug * us%L_to_m**2*g%areaT(i,j) * &
602 (iob%salt_flux(i-i0,j-j0) / cs%ice_salt_concentration)
603 net_fw2(i,j) = net_fw(i,j) / (us%L_to_m**2*g%areaT(i,j))
606 if (cs%adjust_net_fresh_water_by_scaling)
then
607 call adjust_area_mean_to_zero(net_fw2, g, fluxes%netFWGlobalScl)
608 do j=js,je ;
do i=is,ie
609 fluxes%vprec(i,j) = fluxes%vprec(i,j) + us%kg_m2s_to_RZ_T * &
610 (net_fw2(i,j) - net_fw(i,j)/(us%L_to_m**2*g%areaT(i,j))) * g%mask2dT(i,j)
613 fluxes%netFWGlobalAdj =
reproducing_sum(net_fw(:,:), isr, ier, jsr, jer) / cs%area_surf
614 do j=js,je ;
do i=is,ie
615 fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion * fluxes%netFWGlobalAdj ) * g%mask2dT(i,j)
622 if (
associated(fluxes%ustar) .and.
associated(fluxes%ustar_gustless))
then
623 call extract_iob_stresses(iob, index_bounds, time, g, us, cs, ustar=fluxes%ustar, &
624 gustless_ustar=fluxes%ustar_gustless)
625 elseif (
associated(fluxes%ustar))
then
626 call extract_iob_stresses(iob, index_bounds, time, g, us, cs, ustar=fluxes%ustar)
627 elseif (
associated(fluxes%ustar_gustless))
then
628 call extract_iob_stresses(iob, index_bounds, time, g, us, cs, gustless_ustar=fluxes%ustar_gustless)
631 if (coupler_type_initialized(fluxes%tr_fluxes) .and. &
632 coupler_type_initialized(iob%fluxes)) &
633 call coupler_type_copy_data(iob%fluxes, fluxes%tr_fluxes)
635 if (cs%allow_flux_adjustments)
then
637 call apply_flux_adjustments(g, us, cs, time, fluxes)
641 call user_alter_forcing(sfc_state, fluxes, time, g, cs%urf_CS)
643 call cpu_clock_end(id_clock_forcing)
645 end subroutine convert_iob_to_fluxes
650 subroutine convert_iob_to_forces(IOB, forces, index_bounds, Time, G, US, CS, dt_forcing, reset_avg)
652 target,
intent(in) :: iob
655 integer,
dimension(4),
intent(in) :: index_bounds
656 type(time_type),
intent(in) :: time
662 real,
optional,
intent(in) :: dt_forcing
666 logical,
optional,
intent(in) :: reset_avg
669 real,
dimension(SZI_(G),SZJ_(G)) :: &
680 integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq, i0, j0
681 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, isr, ier, jsr, jer
682 integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
684 call cpu_clock_begin(id_clock_forcing)
686 isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
687 jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
688 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
689 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
690 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
691 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
692 isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
693 i0 = is - isc_bnd ; j0 = js - jsc_bnd
697 if (.not.forces%initialized)
then
701 call safe_alloc_ptr(forces%p_surf,isd,ied,jsd,jed)
702 call safe_alloc_ptr(forces%p_surf_full,isd,ied,jsd,jed)
703 if (cs%use_limited_P_SSH)
then
704 forces%p_surf_SSH => forces%p_surf
706 forces%p_surf_SSH => forces%p_surf_full
709 if (cs%rigid_sea_ice)
then
710 call safe_alloc_ptr(forces%rigidity_ice_u,isdb,iedb,jsd,jed)
711 call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,jsdb,jedb)
714 forces%initialized = .true.
717 if ( (
associated(iob%area_berg) .and. (.not.
associated(forces%area_berg))) .or. &
718 (
associated(iob%mass_berg) .and. (.not.
associated(forces%mass_berg))) ) &
721 if (
associated(iob%ice_rigidity))
then
722 rigidity_at_h(:,:) = 0.0
723 call safe_alloc_ptr(forces%rigidity_ice_u,isdb,iedb,jsd,jed)
724 call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,jsdb,jedb)
727 forces%accumulate_rigidity = .true.
728 if (
associated(forces%rigidity_ice_u)) forces%rigidity_ice_u(:,:) = 0.0
729 if (
associated(forces%rigidity_ice_v)) forces%rigidity_ice_v(:,:) = 0.0
732 if (
present(reset_avg))
then ;
if (reset_avg) forces%dt_force_accum = 0.0 ;
endif
733 wt1 = 0.0 ; wt2 = 1.0
734 if (
present(dt_forcing))
then
735 if ((forces%dt_force_accum > 0.0) .and. (dt_forcing > 0.0))
then
736 wt1 = forces%dt_force_accum / (forces%dt_force_accum + dt_forcing)
739 if (dt_forcing > 0.0)
then
740 forces%dt_force_accum = max(forces%dt_force_accum, 0.0) + dt_forcing
742 forces%dt_force_accum = 0.0
745 forces%dt_force_accum = 0.0
749 if (
associated(iob%p))
then
750 if (cs%max_p_surf >= 0.0)
then
751 do j=js,je ;
do i=is,ie
752 forces%p_surf_full(i,j) = g%mask2dT(i,j) * us%kg_m3_to_R*us%m_s_to_L_T**2*iob%p(i-i0,j-j0)
753 forces%p_surf(i,j) = min(forces%p_surf_full(i,j),cs%max_p_surf)
756 do j=js,je ;
do i=is,ie
757 forces%p_surf_full(i,j) = g%mask2dT(i,j) * us%kg_m3_to_R*us%m_s_to_L_T**2*iob%p(i-i0,j-j0)
758 forces%p_surf(i,j) = forces%p_surf_full(i,j)
762 do j=js,je ;
do i=is,ie
763 forces%p_surf_full(i,j) = 0.0
764 forces%p_surf(i,j) = 0.0
767 forces%accumulate_p_surf = .true.
771 call extract_iob_stresses(iob, index_bounds, time, g, us, cs, taux=forces%taux, tauy=forces%tauy, &
772 ustar=forces%ustar, tau_halo=1)
774 call extract_iob_stresses(iob, index_bounds, time, g, us, cs, taux=forces%taux, tauy=forces%tauy, &
775 ustar=ustar_tmp, tau_halo=1)
776 do j=js,je ;
do i=is,ie
777 forces%ustar(i,j) = wt1*forces%ustar(i,j) + wt2*ustar_tmp(i,j)
782 if (cs%approx_net_mass_src .and.
associated(forces%net_mass_src))
then
783 net_mass_src(:,:) = 0.0
784 i0 = is - isc_bnd ; j0 = js - jsc_bnd
785 do j=js,je ;
do i=is,ie ;
if (g%mask2dT(i,j) > 0.0)
then
786 if (
associated(iob%lprec)) &
787 net_mass_src(i,j) = net_mass_src(i,j) + iob%lprec(i-i0,j-j0)
788 if (
associated(iob%fprec)) &
789 net_mass_src(i,j) = net_mass_src(i,j) + iob%fprec(i-i0,j-j0)
790 if (
associated(iob%runoff)) &
791 net_mass_src(i,j) = net_mass_src(i,j) + iob%runoff(i-i0,j-j0)
792 if (
associated(iob%calving)) &
793 net_mass_src(i,j) = net_mass_src(i,j) + iob%calving(i-i0,j-j0)
794 if (
associated(iob%q_flux)) &
795 net_mass_src(i,j) = net_mass_src(i,j) - iob%q_flux(i-i0,j-j0)
796 endif ;
enddo ;
enddo
798 do j=js,je ;
do i=is,ie
799 forces%net_mass_src(i,j) = wt2*net_mass_src(i,j)
802 do j=js,je ;
do i=is,ie
803 forces%net_mass_src(i,j) = wt1*forces%net_mass_src(i,j) + wt2*net_mass_src(i,j)
806 forces%net_mass_src_set = .true.
808 forces%net_mass_src_set = .false.
812 if (
associated(iob%area_berg))
then ;
do j=js,je ;
do i=is,ie
813 forces%area_berg(i,j) = iob%area_berg(i-i0,j-j0) * g%mask2dT(i,j)
814 enddo ;
enddo ;
endif
816 if (
associated(iob%mass_berg))
then ;
do j=js,je ;
do i=is,ie
817 forces%mass_berg(i,j) = us%m_to_Z*us%kg_m3_to_R * iob%mass_berg(i-i0,j-j0) * g%mask2dT(i,j)
818 enddo ;
enddo ;
endif
821 if (
associated(iob%ice_rigidity))
then
822 do j=js,je ;
do i=is,ie
823 rigidity_at_h(i,j) = us%m_to_L**3*us%Z_to_L*us%T_to_s * iob%ice_rigidity(i-i0,j-j0) * g%mask2dT(i,j)
825 call pass_var(rigidity_at_h, g%Domain, halo=1)
826 do i=is-1,ie ;
do j=js,je
827 forces%rigidity_ice_u(i,j) = forces%rigidity_ice_u(i,j) + &
828 min(rigidity_at_h(i,j), rigidity_at_h(i+1,j))
830 do i=is,ie ;
do j=js-1,je
831 forces%rigidity_ice_v(i,j) = forces%rigidity_ice_v(i,j) + &
832 min(rigidity_at_h(i,j), rigidity_at_h(i,j+1))
836 if (cs%rigid_sea_ice)
then
837 call pass_var(forces%p_surf_full, g%Domain, halo=1)
838 i_gearth = 1.0 / cs%g_Earth
839 kv_rho_ice = (cs%Kv_sea_ice / cs%density_sea_ice)
840 do i=is-1,ie ;
do j=js,je
841 mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i+1,j)) * i_gearth
843 if (mass_ice > cs%rigid_sea_ice_mass)
then
844 mass_eff = (mass_ice - cs%rigid_sea_ice_mass)**2 / (mass_ice + cs%rigid_sea_ice_mass)
846 forces%rigidity_ice_u(i,j) = forces%rigidity_ice_u(i,j) + kv_rho_ice * mass_eff
848 do i=is,ie ;
do j=js-1,je
849 mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i,j+1)) * i_gearth
851 if (mass_ice > cs%rigid_sea_ice_mass)
then
852 mass_eff = (mass_ice - cs%rigid_sea_ice_mass)**2 / (mass_ice + cs%rigid_sea_ice_mass)
854 forces%rigidity_ice_v(i,j) = forces%rigidity_ice_v(i,j) + kv_rho_ice * mass_eff
858 if (cs%allow_flux_adjustments)
then
860 call apply_force_adjustments(g, us, cs, time, forces)
866 call cpu_clock_end(id_clock_forcing)
867 end subroutine convert_iob_to_forces
873 subroutine extract_iob_stresses(IOB, index_bounds, Time, G, US, CS, taux, tauy, ustar, &
874 gustless_ustar, tau_halo)
876 target,
intent(in) :: iob
878 integer,
dimension(4),
intent(in) :: index_bounds
879 type(time_type),
intent(in) :: time
885 real,
dimension(SZIB_(G),SZJ_(G)), &
886 optional,
intent(inout) :: taux
887 real,
dimension(SZI_(G),SZJB_(G)), &
888 optional,
intent(inout) :: tauy
889 real,
dimension(SZI_(G),SZJ_(G)), &
890 optional,
intent(inout) :: ustar
891 real,
dimension(SZI_(G),SZJ_(G)), &
892 optional,
intent(out) :: gustless_ustar
894 integer,
optional,
intent(in) :: tau_halo
897 real,
dimension(SZI_(G),SZJ_(G)) :: taux_in_a
898 real,
dimension(SZI_(G),SZJ_(G)) :: tauy_in_a
899 real,
dimension(SZIB_(G),SZJ_(G)) :: taux_in_c
900 real,
dimension(SZI_(G),SZJB_(G)) :: tauy_in_c
901 real,
dimension(SZIB_(G),SZJB_(G)) :: taux_in_b
902 real,
dimension(SZIB_(G),SZJB_(G)) :: tauy_in_b
908 real :: pa_conversion
909 real :: stress_conversion
911 logical :: do_ustar, do_gustless
912 integer :: wind_stagger
913 integer :: i, j, is, ie, js, je, ish, ieh, jsh, jeh, isqh, ieqh, jsqh, jeqh, i0, j0, halo
915 halo = 0 ;
if (
present(tau_halo)) halo = tau_halo
916 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
917 ish = g%isc-halo ; ieh = g%iec+halo ; jsh = g%jsc-halo ; jeh = g%jec+halo
918 isqh = g%IscB-halo ; ieqh = g%IecB+halo ; jsqh = g%JscB-halo ; jeqh = g%JecB+halo
919 i0 = is - index_bounds(1) ; j0 = js - index_bounds(3)
921 irho0 = us%L_to_Z / cs%Rho0
922 pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
923 stress_conversion = pa_conversion * cs%wind_stress_multiplier
925 do_ustar =
present(ustar) ; do_gustless =
present(gustless_ustar)
927 wind_stagger = cs%wind_stagger
928 if ((iob%wind_stagger == agrid) .or. (iob%wind_stagger == bgrid_ne) .or. &
929 (iob%wind_stagger == cgrid_ne)) wind_stagger = iob%wind_stagger
931 if (
associated(iob%u_flux).neqv.
associated(iob%v_flux))
call mom_error(fatal,
"extract_IOB_stresses: "//&
932 "associated(IOB%u_flux) /= associated(IOB%v_flux !!!")
933 if (
present(taux).neqv.
present(tauy))
call mom_error(fatal,
"extract_IOB_stresses: "//&
934 "present(taux) /= present(tauy) !!!")
937 if (
present(taux) .or.
present(tauy) .or. &
938 ((do_ustar.or.do_gustless) .and. .not.
associated(iob%stress_mag)) )
then
940 if (wind_stagger == bgrid_ne)
then
941 taux_in_b(:,:) = 0.0 ; tauy_in_b(:,:) = 0.0
942 if (
associated(iob%u_flux).and.
associated(iob%v_flux))
then
943 do j=js,je ;
do i=is,ie
944 taux_in_b(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
945 tauy_in_b(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
950 call pass_vector(taux_in_b, tauy_in_b, g%Domain, stagger=bgrid_ne, halo=max(1,halo))
952 if (
present(taux).and.
present(tauy))
then
953 do j=jsh,jeh ;
do i=isqh,ieqh
955 if ((g%mask2dBu(i,j) + g%mask2dBu(i,j-1)) > 0) &
956 taux(i,j) = (g%mask2dBu(i,j)*taux_in_b(i,j) + g%mask2dBu(i,j-1)*taux_in_b(i,j-1)) / &
957 (g%mask2dBu(i,j) + g%mask2dBu(i,j-1))
959 do j=jsqh,jeqh ;
do i=ish,ieh
961 if ((g%mask2dBu(i,j) + g%mask2dBu(i-1,j)) > 0) &
962 tauy(i,j) = (g%mask2dBu(i,j)*tauy_in_b(i,j) + g%mask2dBu(i-1,j)*tauy_in_b(i-1,j)) / &
963 (g%mask2dBu(i,j) + g%mask2dBu(i-1,j))
966 elseif (wind_stagger == agrid)
then
967 taux_in_a(:,:) = 0.0 ; tauy_in_a(:,:) = 0.0
968 if (
associated(iob%u_flux).and.
associated(iob%v_flux))
then
969 do j=js,je ;
do i=is,ie
970 taux_in_a(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
971 tauy_in_a(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
976 call pass_vector(taux_in_a, tauy_in_a, g%Domain, to_all+omit_corners, stagger=agrid, halo=1)
978 call pass_vector(taux_in_a, tauy_in_a, g%Domain, stagger=agrid, halo=max(1,halo))
981 if (
present(taux))
then ;
do j=jsh,jeh ;
do i=isqh,ieqh
983 if ((g%mask2dT(i,j) + g%mask2dT(i+1,j)) > 0) &
984 taux(i,j) = (g%mask2dT(i,j)*taux_in_a(i,j) + g%mask2dT(i+1,j)*taux_in_a(i+1,j)) / &
985 (g%mask2dT(i,j) + g%mask2dT(i+1,j))
986 enddo ;
enddo ;
endif
988 if (
present(tauy))
then ;
do j=jsqh,jeqh ;
do i=ish,ieh
990 if ((g%mask2dT(i,j) + g%mask2dT(i,j+1)) > 0) &
991 tauy(i,j) = (g%mask2dT(i,j)*tauy_in_a(i,j) + g%mask2dT(i,j+1)*tauy_in_a(i,j+1)) / &
992 (g%mask2dT(i,j) + g%mask2dT(i,j+1))
993 enddo ;
enddo ;
endif
996 taux_in_c(:,:) = 0.0 ; tauy_in_c(:,:) = 0.0
997 if (
associated(iob%u_flux).and.
associated(iob%v_flux))
then
998 do j=js,je ;
do i=is,ie
999 taux_in_c(i,j) = iob%u_flux(i-i0,j-j0) * stress_conversion
1000 tauy_in_c(i,j) = iob%v_flux(i-i0,j-j0) * stress_conversion
1005 call pass_vector(taux_in_c, tauy_in_c, g%Domain, halo=max(1,halo))
1007 if (
present(taux).and.
present(tauy))
then
1008 do j=jsh,jeh ;
do i=isqh,ieqh
1009 taux(i,j) = g%mask2dCu(i,j)*taux_in_c(i,j)
1011 do j=jsqh,jeqh ;
do i=ish,ieh
1012 tauy(i,j) = g%mask2dCv(i,j)*tauy_in_c(i,j)
1018 if (do_ustar .or. do_gustless)
then
1023 if (
associated(iob%stress_mag))
then
1024 if (do_ustar)
then ;
do j=js,je ;
do i=is,ie
1025 gustiness = cs%gust_const
1026 if (cs%read_gust_2d)
then
1027 if ((wind_stagger == cgrid_ne) .or. &
1028 ((wind_stagger == agrid) .and. (g%mask2dT(i,j) > 0)) .or. &
1029 ((wind_stagger == bgrid_ne) .and. &
1030 (((g%mask2dBu(i,j) + g%mask2dBu(i-1,j-1)) + &
1031 (g%mask2dBu(i,j-1) + g%mask2dBu(i-1,j))) > 0)) ) &
1032 gustiness = cs%gust(i,j)
1034 ustar(i,j) = sqrt(gustiness*irho0 + irho0*pa_conversion*iob%stress_mag(i-i0,j-j0))
1035 enddo ;
enddo ;
endif
1036 if (cs%answers_2018)
then
1037 if (do_gustless)
then ;
do j=js,je ;
do i=is,ie
1038 gustless_ustar(i,j) = sqrt(pa_conversion*us%L_to_Z*iob%stress_mag(i-i0,j-j0) / cs%Rho0)
1039 enddo ;
enddo ;
endif
1041 if (do_gustless)
then ;
do j=js,je ;
do i=is,ie
1042 gustless_ustar(i,j) = sqrt(irho0 * pa_conversion*iob%stress_mag(i-i0,j-j0))
1043 enddo ;
enddo ;
endif
1045 elseif (wind_stagger == bgrid_ne)
then
1046 do j=js,je ;
do i=is,ie
1047 tau_mag = 0.0 ; gustiness = cs%gust_const
1048 if (((g%mask2dBu(i,j) + g%mask2dBu(i-1,j-1)) + &
1049 (g%mask2dBu(i,j-1) + g%mask2dBu(i-1,j))) > 0)
then
1050 tau_mag = sqrt(((g%mask2dBu(i,j)*(taux_in_b(i,j)**2 + tauy_in_b(i,j)**2) + &
1051 g%mask2dBu(i-1,j-1)*(taux_in_b(i-1,j-1)**2 + tauy_in_b(i-1,j-1)**2)) + &
1052 (g%mask2dBu(i,j-1)*(taux_in_b(i,j-1)**2 + tauy_in_b(i,j-1)**2) + &
1053 g%mask2dBu(i-1,j)*(taux_in_b(i-1,j)**2 + tauy_in_b(i-1,j)**2)) ) / &
1054 ((g%mask2dBu(i,j) + g%mask2dBu(i-1,j-1)) + (g%mask2dBu(i,j-1) + g%mask2dBu(i-1,j))) )
1055 if (cs%read_gust_2d) gustiness = cs%gust(i,j)
1057 if (do_ustar) ustar(i,j) = sqrt(gustiness*irho0 + irho0 * tau_mag)
1058 if (cs%answers_2018)
then
1059 if (do_gustless) gustless_ustar(i,j) = sqrt(us%L_to_Z*tau_mag / cs%Rho0)
1061 if (do_gustless) gustless_ustar(i,j) = sqrt(irho0 * tau_mag)
1064 elseif (wind_stagger == agrid)
then
1065 do j=js,je ;
do i=is,ie
1066 tau_mag = g%mask2dT(i,j) * sqrt(taux_in_a(i,j)**2 + tauy_in_a(i,j)**2)
1067 gustiness = cs%gust_const
1068 if (cs%read_gust_2d .and. (g%mask2dT(i,j) > 0)) gustiness = cs%gust(i,j)
1069 if (do_ustar) ustar(i,j) = sqrt(gustiness*irho0 + irho0 * tau_mag)
1070 if (cs%answers_2018)
then
1071 if (do_gustless) gustless_ustar(i,j) = sqrt(us%L_to_Z*tau_mag / cs%Rho0)
1073 if (do_gustless) gustless_ustar(i,j) = sqrt(irho0 * tau_mag)
1077 do j=js,je ;
do i=is,ie
1078 taux2 = 0.0 ; tauy2 = 0.0
1079 if ((g%mask2dCu(i-1,j) + g%mask2dCu(i,j)) > 0) &
1080 taux2 = (g%mask2dCu(i-1,j)*taux_in_c(i-1,j)**2 + g%mask2dCu(i,j)*taux_in_c(i,j)**2) / &
1081 (g%mask2dCu(i-1,j) + g%mask2dCu(i,j))
1082 if ((g%mask2dCv(i,j-1) + g%mask2dCv(i,j)) > 0) &
1083 tauy2 = (g%mask2dCv(i,j-1)*tauy_in_c(i,j-1)**2 + g%mask2dCv(i,j)*tauy_in_c(i,j)**2) / &
1084 (g%mask2dCv(i,j-1) + g%mask2dCv(i,j))
1085 tau_mag = sqrt(taux2 + tauy2)
1087 gustiness = cs%gust_const
1088 if (cs%read_gust_2d) gustiness = cs%gust(i,j)
1090 if (do_ustar) ustar(i,j) = sqrt(gustiness*irho0 + irho0 * tau_mag)
1091 if (cs%answers_2018)
then
1092 if (do_gustless) gustless_ustar(i,j) = sqrt(us%L_to_Z*tau_mag / cs%Rho0)
1094 if (do_gustless) gustless_ustar(i,j) = sqrt(irho0 * tau_mag)
1100 end subroutine extract_iob_stresses
1109 subroutine apply_flux_adjustments(G, US, CS, Time, fluxes)
1113 type(time_type),
intent(in) :: time
1114 type(
forcing),
intent(inout) :: fluxes
1117 real,
dimension(SZI_(G),SZJ_(G)) :: temp_at_h
1119 integer :: isc, iec, jsc, jec, i, j
1120 logical :: overrode_h
1122 isc = g%isc; iec = g%iec ; jsc = g%jsc; jec = g%jec
1124 overrode_h = .false.
1125 call data_override(
'OCN',
'hflx_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
1127 if (overrode_h)
then ;
do j=jsc,jec ;
do i=isc,iec
1128 fluxes%heat_added(i,j) = fluxes%heat_added(i,j) + us%W_m2_to_QRZ_T*temp_at_h(i,j)* g%mask2dT(i,j)
1129 enddo ;
enddo ;
endif
1132 overrode_h = .false.
1133 call data_override(
'OCN',
'sflx_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
1135 if (overrode_h)
then ;
do j=jsc,jec ;
do i=isc,iec
1136 fluxes%salt_flux_added(i,j) = fluxes%salt_flux_added(i,j) + &
1137 us%kg_m2s_to_RZ_T * temp_at_h(i,j)* g%mask2dT(i,j)
1138 enddo ;
enddo ;
endif
1141 overrode_h = .false.
1142 call data_override(
'OCN',
'prcme_adj', temp_at_h(isc:iec,jsc:jec), time, override=overrode_h)
1144 if (overrode_h)
then ;
do j=jsc,jec ;
do i=isc,iec
1145 fluxes%vprec(i,j) = fluxes%vprec(i,j) + us%kg_m2s_to_RZ_T * temp_at_h(i,j)* g%mask2dT(i,j)
1146 enddo ;
enddo ;
endif
1148 end subroutine apply_flux_adjustments
1155 subroutine apply_force_adjustments(G, US, CS, Time, forces)
1159 type(time_type),
intent(in) :: time
1163 real,
dimension(SZI_(G),SZJ_(G)) :: tempx_at_h
1164 real,
dimension(SZI_(G),SZJ_(G)) :: tempy_at_h
1166 integer :: isc, iec, jsc, jec, i, j
1167 real :: dlondx, dlondy, rdlon, cosa, sina, zonal_tau, merid_tau
1168 real :: pa_conversion
1169 logical :: overrode_x, overrode_y
1171 isc = g%isc; iec = g%iec ; jsc = g%jsc; jec = g%jec
1172 pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
1174 tempx_at_h(:,:) = 0.0 ; tempy_at_h(:,:) = 0.0
1176 overrode_x = .false. ; overrode_y = .false.
1177 call data_override(
'OCN',
'taux_adj', tempx_at_h(isc:iec,jsc:jec), time, override=overrode_x)
1178 call data_override(
'OCN',
'tauy_adj', tempy_at_h(isc:iec,jsc:jec), time, override=overrode_y)
1180 if (overrode_x .or. overrode_y)
then
1181 if (.not. (overrode_x .and. overrode_y))
call mom_error(fatal,
"apply_flux_adjustments: "//&
1182 "Both taux_adj and tauy_adj must be specified, or neither, in data_table")
1185 call pass_vector(tempx_at_h, tempy_at_h, g%Domain, to_all, agrid, halo=1)
1186 do j=jsc-1,jec+1 ;
do i=isc-1,iec+1
1187 dlondx = g%geoLonCu(i,j) - g%geoLonCu(i-1,j)
1188 dlondy = g%geoLonCv(i,j) - g%geoLonCv(i,j-1)
1189 rdlon = sqrt( dlondx * dlondx + dlondy * dlondy )
1190 if (rdlon > 0.) rdlon = 1. / rdlon
1191 cosa = dlondx * rdlon
1192 sina = dlondy * rdlon
1193 zonal_tau = pa_conversion * tempx_at_h(i,j)
1194 merid_tau = pa_conversion * tempy_at_h(i,j)
1195 tempx_at_h(i,j) = cosa * zonal_tau - sina * merid_tau
1196 tempy_at_h(i,j) = sina * zonal_tau + cosa * merid_tau
1200 do j=jsc,jec ;
do i=isc-1,iec
1201 forces%taux(i,j) = forces%taux(i,j) + 0.5 * ( tempx_at_h(i,j) + tempx_at_h(i+1,j) )
1204 do j=jsc-1,jec ;
do i=isc,iec
1205 forces%tauy(i,j) = forces%tauy(i,j) + 0.5 * ( tempy_at_h(i,j) + tempy_at_h(i,j+1) )
1209 end subroutine apply_force_adjustments
1212 subroutine forcing_save_restart(CS, G, Time, directory, time_stamped, &
1217 type(time_type),
intent(in) :: time
1218 character(len=*),
intent(in) :: directory
1220 logical,
optional,
intent(in) :: time_stamped
1222 character(len=*),
optional,
intent(in) :: filename_suffix
1225 if (.not.
associated(cs))
return
1226 if (.not.
associated(cs%restart_CSp))
return
1227 call save_restart(directory, time, g, cs%restart_CSp, time_stamped)
1229 end subroutine forcing_save_restart
1232 subroutine surface_forcing_init(Time, G, US, param_file, diag, CS, wind_stagger)
1233 type(time_type),
intent(in) :: time
1237 type(
diag_ctrl),
target,
intent(inout) :: diag
1241 integer,
optional,
intent(in) :: wind_stagger
1247 logical :: new_sim, iceberg_flux_diags
1248 logical :: default_2018_answers
1249 type(time_type) :: time_frc
1250 character(len=200) :: tideamp_file, gust_file, salt_file, temp_file
1252 # include "version_variable.h"
1253 character(len=40) :: mdl =
"MOM_surface_forcing"
1254 character(len=48) :: stagger
1255 character(len=48) :: flnam
1256 character(len=240) :: basin_file
1257 integer :: i, j, isd, ied, jsd, jed
1259 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1261 if (
associated(cs))
then
1262 call mom_error(warning,
"surface_forcing_init called with an associated "// &
1263 "control structure.")
1268 id_clock_forcing=cpu_clock_id(
'Ocean surface forcing', grain=clock_subcomponent)
1269 call cpu_clock_begin(id_clock_forcing)
1273 call write_version_number(version)
1275 call log_version(param_file, mdl, version,
"", log_to_all=.true., debugging=.true.)
1277 call get_param(param_file, mdl,
"INPUTDIR", cs%inputdir, &
1278 "The directory in which all input files are found.", &
1280 cs%inputdir = slasher(cs%inputdir)
1281 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", cs%use_temperature, &
1282 "If true, Temperature and salinity are used as state "//&
1283 "variables.", default=.true.)
1284 call get_param(param_file, mdl,
"RHO_0", cs%Rho0, &
1285 "The mean ocean density used with BOUSSINESQ true to "//&
1286 "calculate accelerations and the mass for conservation "//&
1287 "properties, or with BOUSSINSEQ false to convert some "//&
1288 "parameters from vertical units of m to kg m-2.", &
1289 units=
"kg m-3", default=1035.0, scale=us%kg_m3_to_R)
1290 call get_param(param_file, mdl,
"LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
1291 "The latent heat of fusion.", units=
"J/kg", default=hlf)
1292 call get_param(param_file, mdl,
"LATENT_HEAT_VAPORIZATION", cs%latent_heat_vapor, &
1293 "The latent heat of fusion.", units=
"J/kg", default=hlv)
1294 call get_param(param_file, mdl,
"MAX_P_SURF", cs%max_p_surf, &
1295 "The maximum surface pressure that can be exerted by the "//&
1296 "atmosphere and floating sea-ice or ice shelves. This is "//&
1297 "needed because the FMS coupling structure does not "//&
1298 "limit the water that can be frozen out of the ocean and "//&
1299 "the ice-ocean heat fluxes are treated explicitly. No "//&
1300 "limit is applied if a negative value is used.", &
1301 units=
"Pa", default=-1.0, scale=us%kg_m3_to_R*us%m_s_to_L_T**2)
1302 call get_param(param_file, mdl,
"RESTORE_SALINITY", cs%restore_salt, &
1303 "If true, the coupled driver will add a globally-balanced "//&
1304 "fresh-water flux that drives sea-surface salinity "//&
1305 "toward specified values.", default=.false.)
1306 call get_param(param_file, mdl,
"RESTORE_TEMPERATURE", cs%restore_temp, &
1307 "If true, the coupled driver will add a "//&
1308 "heat flux that drives sea-surface temperature "//&
1309 "toward specified values.", default=.false.)
1310 call get_param(param_file, mdl,
"ADJUST_NET_SRESTORE_TO_ZERO", &
1311 cs%adjust_net_srestore_to_zero, &
1312 "If true, adjusts the salinity restoring seen to zero "//&
1313 "whether restoring is via a salt flux or virtual precip.",&
1314 default=cs%restore_salt)
1315 call get_param(param_file, mdl,
"ADJUST_NET_SRESTORE_BY_SCALING", &
1316 cs%adjust_net_srestore_by_scaling, &
1317 "If true, adjustments to salt restoring to achieve zero net are "//&
1318 "made by scaling values without moving the zero contour.",&
1320 call get_param(param_file, mdl,
"ADJUST_NET_FRESH_WATER_TO_ZERO", &
1321 cs%adjust_net_fresh_water_to_zero, &
1322 "If true, adjusts the net fresh-water forcing seen "//&
1323 "by the ocean (including restoring) to zero.", default=.false.)
1324 if (cs%adjust_net_fresh_water_to_zero) &
1325 call get_param(param_file, mdl,
"USE_NET_FW_ADJUSTMENT_SIGN_BUG", &
1326 cs%use_net_FW_adjustment_sign_bug, &
1327 "If true, use the wrong sign for the adjustment to "//&
1328 "the net fresh-water.", default=.false.)
1329 call get_param(param_file, mdl,
"ADJUST_NET_FRESH_WATER_BY_SCALING", &
1330 cs%adjust_net_fresh_water_by_scaling, &
1331 "If true, adjustments to net fresh water to achieve zero net are "//&
1332 "made by scaling values without moving the zero contour.",&
1334 call get_param(param_file, mdl,
"ICE_SALT_CONCENTRATION", &
1335 cs%ice_salt_concentration, &
1336 "The assumed sea-ice salinity needed to reverse engineer the "//&
1337 "melt flux (or ice-ocean fresh-water flux).", &
1338 units=
"kg/kg", default=0.005)
1339 call get_param(param_file, mdl,
"USE_LIMITED_PATM_SSH", cs%use_limited_P_SSH, &
1340 "If true, return the sea surface height with the "//&
1341 "correction for the atmospheric (and sea-ice) pressure "//&
1342 "limited by max_p_surf instead of the full atmospheric "//&
1343 "pressure.", default=.true.)
1344 call get_param(param_file, mdl,
"APPROX_NET_MASS_SRC", cs%approx_net_mass_src, &
1345 "If true, use the net mass sources from the ice-ocean "//&
1346 "boundary type without any further adjustments to drive "//&
1347 "the ocean dynamics. The actual net mass source may differ "//&
1348 "due to internal corrections.", default=.false.)
1350 if (
present(wind_stagger))
then
1351 if (wind_stagger == agrid)
then ; stagger =
'AGRID'
1352 elseif (wind_stagger == bgrid_ne)
then ; stagger =
'BGRID_NE'
1353 elseif (wind_stagger == cgrid_ne)
then ; stagger =
'CGRID_NE'
1354 else ; stagger =
'UNKNOWN' ;
call mom_error(fatal,
"surface_forcing_init: WIND_STAGGER = "// &
1355 trim(stagger)//
"is invalid.");
endif
1356 call log_param(param_file, mdl,
"WIND_STAGGER", stagger, &
1357 "The staggering of the input wind stress field "//&
1358 "from the coupler that is actually used.")
1359 cs%wind_stagger = wind_stagger
1361 call get_param(param_file, mdl,
"WIND_STAGGER", stagger, &
1362 "A case-insensitive character string to indicate the "//&
1363 "staggering of the input wind stress field. Valid "//&
1364 "values are 'A', 'B', or 'C'.", default=
"C")
1365 if (uppercase(stagger(1:1)) ==
'A')
then ; cs%wind_stagger = agrid
1366 elseif (uppercase(stagger(1:1)) ==
'B')
then ; cs%wind_stagger = bgrid_ne
1367 elseif (uppercase(stagger(1:1)) ==
'C')
then ; cs%wind_stagger = cgrid_ne
1368 else ;
call mom_error(fatal,
"surface_forcing_init: WIND_STAGGER = "// &
1369 trim(stagger)//
" is invalid.") ;
endif
1372 call get_param(param_file, mdl,
"WIND_STRESS_MULTIPLIER", cs%wind_stress_multiplier, &
1373 "A factor multiplying the wind-stress given to the ocean by the "//&
1374 "coupler. This is used for testing and should be =1.0 for any "//&
1375 "production runs.", units=
"nondim", default=1.0)
1377 if (cs%restore_salt)
then
1378 call get_param(param_file, mdl,
"FLUXCONST", cs%Flux_const, &
1379 "The constant that relates the restoring surface fluxes to the relative "//&
1380 "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
1381 default=0.0, units=
"m day-1", scale=us%m_to_Z*us%T_to_s)
1383 cs%Flux_const = cs%Flux_const / 86400.0
1384 call get_param(param_file, mdl,
"SALT_RESTORE_FILE", cs%salt_restore_file, &
1385 "A file in which to find the surface salinity to use for restoring.", &
1386 default=
"salt_restore.nc")
1387 call get_param(param_file, mdl,
"SALT_RESTORE_VARIABLE", cs%salt_restore_var_name, &
1388 "The name of the surface salinity variable to read from "//&
1389 "SALT_RESTORE_FILE for restoring salinity.", &
1392 call get_param(param_file, mdl,
"SRESTORE_AS_SFLUX", cs%salt_restore_as_sflux, &
1393 "If true, the restoring of salinity is applied as a salt "//&
1394 "flux instead of as a freshwater flux.", default=.false.)
1395 call get_param(param_file, mdl,
"MAX_DELTA_SRESTORE", cs%max_delta_srestore, &
1396 "The maximum salinity difference used in restoring terms.", &
1397 units=
"PSU or g kg-1", default=999.0)
1398 call get_param(param_file, mdl,
"MASK_SRESTORE_UNDER_ICE", &
1399 cs%mask_srestore_under_ice, &
1400 "If true, disables SSS restoring under sea-ice based on a frazil "//&
1401 "criteria (SST<=Tf). Only used when RESTORE_SALINITY is True.", &
1403 call get_param(param_file, mdl,
"MASK_SRESTORE_MARGINAL_SEAS", &
1404 cs%mask_srestore_marginal_seas, &
1405 "If true, disable SSS restoring in marginal seas. Only used when "//&
1406 "RESTORE_SALINITY is True.", default=.false.)
1407 call get_param(param_file, mdl,
"BASIN_FILE", basin_file, &
1408 "A file in which to find the basin masks, in variable 'basin'.", &
1410 basin_file = trim(cs%inputdir) // trim(basin_file)
1411 call safe_alloc_ptr(cs%basin_mask,isd,ied,jsd,jed) ; cs%basin_mask(:,:) = 1.0
1412 if (cs%mask_srestore_marginal_seas)
then
1413 call mom_read_data(basin_file,
'basin',cs%basin_mask,g%domain, timelevel=1)
1414 do j=jsd,jed ;
do i=isd,ied
1415 if (cs%basin_mask(i,j) >= 6.0)
then ; cs%basin_mask(i,j) = 0.0
1416 else ; cs%basin_mask(i,j) = 1.0 ;
endif
1419 call get_param(param_file, mdl,
"MASK_SRESTORE", cs%mask_srestore, &
1420 "If true, read a file (salt_restore_mask) containing "//&
1421 "a mask for SSS restoring.", default=.false.)
1424 if (cs%restore_temp)
then
1425 call get_param(param_file, mdl,
"FLUXCONST", cs%Flux_const, &
1426 "The constant that relates the restoring surface fluxes to the relative "//&
1427 "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
1428 default=0.0, units=
"m day-1", scale=us%m_to_Z*us%T_to_s)
1430 cs%Flux_const = cs%Flux_const / 86400.0
1431 call get_param(param_file, mdl,
"SST_RESTORE_FILE", cs%temp_restore_file, &
1432 "A file in which to find the surface temperature to use for restoring.", &
1433 default=
"temp_restore.nc")
1434 call get_param(param_file, mdl,
"SST_RESTORE_VARIABLE", cs%temp_restore_var_name, &
1435 "The name of the surface temperature variable to read from "//&
1436 "SST_RESTORE_FILE for restoring sst.", &
1439 call get_param(param_file, mdl,
"MAX_DELTA_TRESTORE", cs%max_delta_trestore, &
1440 "The maximum sst difference used in restoring terms.", &
1441 units=
"degC ", default=999.0)
1442 call get_param(param_file, mdl,
"MASK_TRESTORE", cs%mask_trestore, &
1443 "If true, read a file (temp_restore_mask) containing "//&
1444 "a mask for SST restoring.", default=.false.)
1452 call get_param(param_file, mdl,
"CD_TIDES", cs%cd_tides, &
1453 "The drag coefficient that applies to the tides.", &
1454 units=
"nondim", default=1.0e-4)
1455 call get_param(param_file, mdl,
"READ_TIDEAMP", cs%read_TIDEAMP, &
1456 "If true, read a file (given by TIDEAMP_FILE) containing "//&
1457 "the tidal amplitude with INT_TIDE_DISSIPATION.", default=.false.)
1458 if (cs%read_TIDEAMP)
then
1459 call get_param(param_file, mdl,
"TIDEAMP_FILE", tideamp_file, &
1460 "The path to the file containing the spatially varying "//&
1461 "tidal amplitudes with INT_TIDE_DISSIPATION.", &
1462 default=
"tideamp.nc")
1465 call get_param(param_file, mdl,
"UTIDE", cs%utide, &
1466 "The constant tidal amplitude used with INT_TIDE_DISSIPATION.", &
1467 units=
"m s-1", default=0.0, scale=us%m_to_Z*us%T_to_s)
1470 call safe_alloc_ptr(cs%TKE_tidal,isd,ied,jsd,jed)
1471 call safe_alloc_ptr(cs%ustar_tidal,isd,ied,jsd,jed)
1473 if (cs%read_TIDEAMP)
then
1474 tideamp_file = trim(cs%inputdir) // trim(tideamp_file)
1475 call mom_read_data(tideamp_file,
'tideamp',cs%TKE_tidal,g%domain,timelevel=1, scale=us%m_to_Z*us%T_to_s)
1476 do j=jsd, jed;
do i=isd, ied
1477 utide = cs%TKE_tidal(i,j)
1478 cs%TKE_tidal(i,j) = g%mask2dT(i,j)*cs%Rho0*cs%cd_tides*(utide*utide*utide)
1479 cs%ustar_tidal(i,j) = sqrt(cs%cd_tides)*utide
1482 do j=jsd,jed;
do i=isd,ied
1484 cs%TKE_tidal(i,j) = cs%Rho0*cs%cd_tides*(utide*utide*utide)
1485 cs%ustar_tidal(i,j) = sqrt(cs%cd_tides)*utide
1489 call time_interp_external_init
1492 call get_param(param_file, mdl,
"READ_GUST_2D", cs%read_gust_2d, &
1493 "If true, use a 2-dimensional gustiness supplied from "//&
1494 "an input file", default=.false.)
1495 call get_param(param_file, mdl,
"GUST_CONST", cs%gust_const, &
1496 "The background gustiness in the winds.", &
1497 units=
"Pa", default=0.0, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1498 if (cs%read_gust_2d)
then
1499 call get_param(param_file, mdl,
"GUST_2D_FILE", gust_file, &
1500 "The file in which the wind gustiness is found in "//&
1501 "variable gustiness.")
1503 call safe_alloc_ptr(cs%gust,isd,ied,jsd,jed)
1504 gust_file = trim(cs%inputdir) // trim(gust_file)
1505 call mom_read_data(gust_file,
'gustiness', cs%gust, g%domain, timelevel=1, &
1506 scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1508 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1509 "This sets the default value for the various _2018_ANSWERS parameters.", &
1511 call get_param(param_file, mdl,
"SURFACE_FORCING_2018_ANSWERS", cs%answers_2018, &
1512 "If true, use the order of arithmetic and expressions that recover the answers "//&
1513 "from the end of 2018. Otherwise, use a simpler expression to calculate gustiness.", &
1514 default=default_2018_answers)
1515 call get_param(param_file, mdl,
"FIX_USTAR_GUSTLESS_BUG", cs%fix_ustar_gustless_bug, &
1516 "If true correct a bug in the time-averaging of the gustless wind friction velocity", &
1520 call get_param(param_file, mdl,
"USE_RIGID_SEA_ICE", cs%rigid_sea_ice, &
1521 "If true, sea-ice is rigid enough to exert a "//&
1522 "nonhydrostatic pressure that resist vertical motion.", &
1524 if (cs%rigid_sea_ice)
then
1525 call get_param(param_file, mdl,
"G_EARTH", cs%g_Earth, &
1526 "The gravitational acceleration of the Earth.", &
1527 units=
"m s-2", default = 9.80, scale=us%Z_to_m*us%m_s_to_L_T**2)
1528 call get_param(param_file, mdl,
"SEA_ICE_MEAN_DENSITY", cs%density_sea_ice, &
1529 "A typical density of sea ice, used with the kinematic "//&
1530 "viscosity, when USE_RIGID_SEA_ICE is true.", &
1531 units=
"kg m-3", default=900.0, scale=us%kg_m3_to_R)
1532 call get_param(param_file, mdl,
"SEA_ICE_VISCOSITY", cs%Kv_sea_ice, &
1533 "The kinematic viscosity of sufficiently thick sea ice "//&
1534 "for use in calculating the rigidity of sea ice.", &
1535 units=
"m2 s-1", default=1.0e9, scale=us%Z_to_L**2*us%m_to_L**2*us%T_to_s)
1536 call get_param(param_file, mdl,
"SEA_ICE_RIGID_MASS", cs%rigid_sea_ice_mass, &
1537 "The mass of sea-ice per unit area at which the sea-ice "//&
1538 "starts to exhibit rigidity", &
1539 units=
"kg m-2", default=1000.0, scale=us%kg_m3_to_R*us%m_to_Z)
1542 call get_param(param_file, mdl,
"ALLOW_ICEBERG_FLUX_DIAGNOSTICS", iceberg_flux_diags, &
1543 "If true, makes available diagnostics of fluxes from icebergs "//&
1544 "as seen by MOM6.", default=.false.)
1545 call register_forcing_type_diags(time, diag, us, cs%use_temperature, cs%handles, &
1546 use_berg_fluxes=iceberg_flux_diags)
1548 call get_param(param_file, mdl,
"ALLOW_FLUX_ADJUSTMENTS", cs%allow_flux_adjustments, &
1549 "If true, allows flux adjustments to specified via the "//&
1550 "data_table using the component name 'OCN'.", default=.false.)
1552 call get_param(param_file, mdl,
"CHECK_NO_LAND_FLUXES", cs%check_no_land_fluxes, &
1553 "If true, checks that values from IOB fluxes are zero "//&
1554 "above land points (i.e. G%mask2dT = 0).", default=.false., &
1555 debuggingparam=.true.)
1557 call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
1559 if (cs%restore_salt)
then
1560 salt_file = trim(cs%inputdir) // trim(cs%salt_restore_file)
1561 cs%id_srestore = init_external_field(salt_file, cs%salt_restore_var_name, domain=g%Domain%mpp_domain)
1562 call safe_alloc_ptr(cs%srestore_mask,isd,ied,jsd,jed); cs%srestore_mask(:,:) = 1.0
1563 if (cs%mask_srestore)
then
1564 flnam = trim(cs%inputdir) //
'salt_restore_mask.nc'
1565 call mom_read_data(flnam,
'mask', cs%srestore_mask, g%domain, timelevel=1)
1569 if (cs%restore_temp)
then
1570 temp_file = trim(cs%inputdir) // trim(cs%temp_restore_file)
1571 cs%id_trestore = init_external_field(temp_file, cs%temp_restore_var_name, domain=g%Domain%mpp_domain)
1572 call safe_alloc_ptr(cs%trestore_mask,isd,ied,jsd,jed); cs%trestore_mask(:,:) = 1.0
1573 if (cs%mask_trestore)
then
1574 flnam = trim(cs%inputdir) //
'temp_restore_mask.nc'
1575 call mom_read_data(flnam,
'mask', cs%trestore_mask, g%domain, timelevel=1)
1580 call restart_init(param_file, cs%restart_CSp,
"MOM_forcing.res")
1583 call restart_init_end(cs%restart_CSp)
1585 if (
associated(cs%restart_CSp))
then
1586 call get_mom_input(dirs=dirs)
1589 if ((dirs%input_filename(1:1) ==
'n') .and. &
1590 (len_trim(dirs%input_filename) == 1)) new_sim = .true.
1591 if (.not.new_sim)
then
1592 call restore_state(dirs%input_filename, dirs%restart_input_dir, time_frc, &
1599 call user_revise_forcing_init(param_file, cs%urf_CS)
1601 call cpu_clock_end(id_clock_forcing)
1602 end subroutine surface_forcing_init
1605 subroutine surface_forcing_end(CS, fluxes)
1609 type(
forcing),
optional,
intent(inout) :: fluxes
1613 if (
present(fluxes))
call deallocate_forcing_type(fluxes)
1617 if (
associated(cs))
deallocate(cs)
1620 end subroutine surface_forcing_end
1623 subroutine ice_ocn_bnd_type_chksum(id, timestep, iobt)
1625 character(len=*),
intent(in) :: id
1626 integer,
intent(in) :: timestep
1630 integer :: n,m, outunit
1634 write(outunit,*)
"BEGIN CHECKSUM(ice_ocean_boundary_type):: ", id, timestep
1635 write(outunit,100)
'iobt%u_flux ', mpp_chksum( iobt%u_flux )
1636 write(outunit,100)
'iobt%v_flux ', mpp_chksum( iobt%v_flux )
1637 write(outunit,100)
'iobt%t_flux ', mpp_chksum( iobt%t_flux )
1638 write(outunit,100)
'iobt%q_flux ', mpp_chksum( iobt%q_flux )
1639 write(outunit,100)
'iobt%salt_flux ', mpp_chksum( iobt%salt_flux )
1640 write(outunit,100)
'iobt%lw_flux ', mpp_chksum( iobt%lw_flux )
1641 write(outunit,100)
'iobt%sw_flux_vis_dir', mpp_chksum( iobt%sw_flux_vis_dir)
1642 write(outunit,100)
'iobt%sw_flux_vis_dif', mpp_chksum( iobt%sw_flux_vis_dif)
1643 write(outunit,100)
'iobt%sw_flux_nir_dir', mpp_chksum( iobt%sw_flux_nir_dir)
1644 write(outunit,100)
'iobt%sw_flux_nir_dif', mpp_chksum( iobt%sw_flux_nir_dif)
1645 write(outunit,100)
'iobt%lprec ', mpp_chksum( iobt%lprec )
1646 write(outunit,100)
'iobt%fprec ', mpp_chksum( iobt%fprec )
1647 write(outunit,100)
'iobt%runoff ', mpp_chksum( iobt%runoff )
1648 write(outunit,100)
'iobt%calving ', mpp_chksum( iobt%calving )
1649 write(outunit,100)
'iobt%p ', mpp_chksum( iobt%p )
1650 if (
associated(iobt%ustar_berg)) &
1651 write(outunit,100)
'iobt%ustar_berg ', mpp_chksum( iobt%ustar_berg )
1652 if (
associated(iobt%area_berg)) &
1653 write(outunit,100)
'iobt%area_berg ', mpp_chksum( iobt%area_berg )
1654 if (
associated(iobt%mass_berg)) &
1655 write(outunit,100)
'iobt%mass_berg ', mpp_chksum( iobt%mass_berg )
1656 100
FORMAT(
" CHECKSUM::",a20,
" = ",z20)
1658 call coupler_type_write_chksums(iobt%fluxes, outunit,
'iobt%')
1660 end subroutine ice_ocn_bnd_type_chksum
1663 subroutine check_mask_val_consistency(val, mask, i, j, varname, G)
1665 real,
intent(in) :: val
1666 real,
intent(in) :: mask
1667 integer,
intent(in) :: i
1668 integer,
intent(in) :: j
1669 character(len=*),
intent(in) :: varname
1672 character(len=48) :: ci, cj
1673 character(len=48) :: ciglo, cjglo
1674 character(len=48) :: cval
1675 character(len=256) :: error_message
1677 if ((mask == 0.) .and. (val /= 0.))
then
1680 write(ciglo,
'(I8)') i + g%HI%idg_offset
1681 write(cjglo,
'(I8)') j + g%HI%jdg_offset
1682 write(cval,
'(E22.16)') val
1683 error_message =
"MOM_surface_forcing: found non-zero value (="//trim(cval)//
") over land "//&
1684 "for variable "//trim(varname)//
" at local point (i, j) = ("//trim(ci)//
", "//trim(cj)//&
1685 ", global point (iglo, jglo) = ("//trim(ciglo)//
", "//trim(cjglo)//
")"
1686 call mom_error(warning, error_message)
1691 end module mom_surface_forcing_gfdl