13 use coupler_types_mod,
only : coupler_1d_bc_type, coupler_2d_bc_type
14 use coupler_types_mod,
only : coupler_type_spawn, coupler_type_destructor
15 use coupler_types_mod,
only : coupler_type_initialized
17 implicit none ;
private
19 #include <MOM_memory.h>
21 public allocate_surface_state, deallocate_surface_state, mom_thermovar_chksum
23 public rotate_surface_state
32 real,
dimension(:,:,:),
pointer :: p => null()
36 real,
dimension(:,:),
pointer :: p => null()
42 real,
allocatable,
dimension(:,:) :: &
43 sst, & !< The sea surface temperature [degC].
44 sss, & !< The sea surface salinity [ppt ~> psu or gSalt/kg].
45 sfc_density, & !< The mixed layer density [R ~> kg m-3].
46 hml, & !< The mixed layer depth [Z ~> m].
47 u, & !< The mixed layer zonal velocity [L T-1 ~> m s-1].
48 v, & !< The mixed layer meridional velocity [L T-1 ~> m s-1].
49 sea_lev, & !< The sea level [Z ~> m]. If a reduced surface gravity is
66 logical :: t_is_cont = .false.
68 logical :: s_is_abss = .false.
70 type(coupler_2d_bc_type) :: tr_fields
74 logical :: arrays_allocated = .false.
82 real,
pointer :: t(:,:,:) => null()
83 real,
pointer :: s(:,:,:) => null()
84 real,
pointer :: p_surf(:,:) => null()
86 type(
eos_type),
pointer :: eqn_of_state => null()
94 logical :: t_is_cont = .false.
96 logical :: s_is_abss = .false.
98 real :: min_salinity = 0.01
101 real,
dimension(:,:),
pointer :: frazil => null()
105 real,
dimension(:,:),
pointer :: salt_deficit => null()
109 real,
dimension(:,:),
pointer :: tempxpme => null()
115 real,
dimension(:,:),
pointer :: internal_heat => null()
121 real,
pointer :: vart(:,:,:) => null()
122 real,
pointer :: vars(:,:,:) => null()
123 real,
pointer :: covarts(:,:,:) => null()
133 real,
pointer,
dimension(:,:,:) :: &
139 real,
pointer,
dimension(:,:,:) :: &
142 real,
pointer,
dimension(:,:,:) :: &
151 u_accel_bt => null(), &
153 real,
pointer,
dimension(:,:,:) :: &
164 real,
pointer,
dimension(:,:,:) :: &
171 du_dt_visc => null(), &
172 dv_dt_visc => null(), &
173 du_dt_dia => null(), &
175 real,
pointer,
dimension(:,:,:) :: du_other => null()
178 real,
pointer,
dimension(:,:,:) :: dv_other => null()
183 real,
pointer :: gradkeu(:,:,:) => null()
184 real,
pointer :: gradkev(:,:,:) => null()
185 real,
pointer :: rv_x_v(:,:,:) => null()
186 real,
pointer :: rv_x_u(:,:,:) => null()
188 real,
pointer :: diag_hfrac_u(:,:,:) => null()
189 real,
pointer :: diag_hfrac_v(:,:,:) => null()
196 real,
pointer,
dimension(:,:,:) :: &
203 real,
pointer :: diapyc_vel(:,:,:) => null()
211 real,
pointer,
dimension(:,:) :: &
212 bbl_thick_u => null(), &
213 bbl_thick_v => null(), &
214 kv_bbl_u => null(), &
215 kv_bbl_v => null(), &
217 real,
pointer,
dimension(:,:) :: tke_bbl => null()
221 real,
pointer,
dimension(:,:) :: &
222 taux_shelf => null(), &
224 real,
pointer,
dimension(:,:) :: tbl_thick_shelf_u => null()
226 real,
pointer,
dimension(:,:) :: tbl_thick_shelf_v => null()
228 real,
pointer,
dimension(:,:) :: kv_tbl_shelf_u => null()
230 real,
pointer,
dimension(:,:) :: kv_tbl_shelf_v => null()
232 real,
pointer,
dimension(:,:) :: nkml_visc_u => null()
237 real,
pointer,
dimension(:,:) :: nkml_visc_v => null()
239 real,
pointer,
dimension(:,:) :: &
241 real,
pointer,
dimension(:,:,:) :: &
244 real,
pointer,
dimension(:,:,:) :: kd_shear => null()
247 real,
pointer,
dimension(:,:,:) :: kv_shear => null()
250 real,
pointer,
dimension(:,:,:) :: kv_shear_bu => null()
253 real,
pointer,
dimension(:,:,:) :: kv_slow => null()
256 real,
pointer,
dimension(:,:,:) :: tke_turb => null()
264 real,
allocatable :: fa_u_ee(:,:)
266 real,
allocatable :: fa_u_e0(:,:)
268 real,
allocatable :: fa_u_w0(:,:)
270 real,
allocatable :: fa_u_ww(:,:)
272 real,
allocatable :: ubt_ww(:,:)
274 real,
allocatable :: ubt_ee(:,:)
276 real,
allocatable :: fa_v_nn(:,:)
278 real,
allocatable :: fa_v_n0(:,:)
280 real,
allocatable :: fa_v_s0(:,:)
282 real,
allocatable :: fa_v_ss(:,:)
284 real,
allocatable :: vbt_ss(:,:)
286 real,
allocatable :: vbt_nn(:,:)
288 real,
allocatable :: h_u(:,:,:)
289 real,
allocatable :: h_v(:,:,:)
290 type(group_pass_type) :: pass_polarity_bt
291 type(group_pass_type) :: pass_fa_uv
298 subroutine allocate_surface_state(sfc_state, G, use_temperature, do_integrals, &
299 gas_fields_ocn, use_meltpot, use_iceshelves, omit_frazil)
301 type(
surface),
intent(inout) :: sfc_state
302 logical,
optional,
intent(in) :: use_temperature
303 logical,
optional,
intent(in) :: do_integrals
305 type(coupler_1d_bc_type), &
306 optional,
intent(in) :: gas_fields_ocn
311 logical,
optional,
intent(in) :: use_meltpot
312 logical,
optional,
intent(in) :: use_iceshelves
314 logical,
optional,
intent(in) :: omit_frazil
318 logical :: use_temp, alloc_integ, use_melt_potential, alloc_iceshelves, alloc_frazil
319 integer :: is, ie, js, je, isd, ied, jsd, jed
320 integer :: isdb, iedb, jsdb, jedb
322 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
323 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
324 isdb = g%isdB ; iedb = g%iedB; jsdb = g%jsdB ; jedb = g%jedB
326 use_temp = .true. ;
if (
present(use_temperature)) use_temp = use_temperature
327 alloc_integ = .true. ;
if (
present(do_integrals)) alloc_integ = do_integrals
328 use_melt_potential = .false. ;
if (
present(use_meltpot)) use_melt_potential = use_meltpot
329 alloc_iceshelves = .false. ;
if (
present(use_iceshelves)) alloc_iceshelves = use_iceshelves
330 alloc_frazil = .true. ;
if (
present(omit_frazil)) alloc_frazil = .not.omit_frazil
332 if (sfc_state%arrays_allocated)
return
335 allocate(sfc_state%SST(isd:ied,jsd:jed)) ; sfc_state%SST(:,:) = 0.0
336 allocate(sfc_state%SSS(isd:ied,jsd:jed)) ; sfc_state%SSS(:,:) = 0.0
338 allocate(sfc_state%sfc_density(isd:ied,jsd:jed)) ; sfc_state%sfc_density(:,:) = 0.0
340 if (use_temp .and. alloc_frazil)
then
341 allocate(sfc_state%frazil(isd:ied,jsd:jed)) ; sfc_state%frazil(:,:) = 0.0
343 allocate(sfc_state%sea_lev(isd:ied,jsd:jed)) ; sfc_state%sea_lev(:,:) = 0.0
344 allocate(sfc_state%Hml(isd:ied,jsd:jed)) ; sfc_state%Hml(:,:) = 0.0
345 allocate(sfc_state%u(isdb:iedb,jsd:jed)) ; sfc_state%u(:,:) = 0.0
346 allocate(sfc_state%v(isd:ied,jsdb:jedb)) ; sfc_state%v(:,:) = 0.0
348 if (use_melt_potential)
then
349 allocate(sfc_state%melt_potential(isd:ied,jsd:jed)) ; sfc_state%melt_potential(:,:) = 0.0
352 if (alloc_integ)
then
354 allocate(sfc_state%ocean_mass(isd:ied,jsd:jed)) ; sfc_state%ocean_mass(:,:) = 0.0
356 allocate(sfc_state%ocean_heat(isd:ied,jsd:jed)) ; sfc_state%ocean_heat(:,:) = 0.0
357 allocate(sfc_state%ocean_salt(isd:ied,jsd:jed)) ; sfc_state%ocean_salt(:,:) = 0.0
358 allocate(sfc_state%TempxPmE(isd:ied,jsd:jed)) ; sfc_state%TempxPmE(:,:) = 0.0
359 allocate(sfc_state%salt_deficit(isd:ied,jsd:jed)) ; sfc_state%salt_deficit(:,:) = 0.0
360 allocate(sfc_state%internal_heat(isd:ied,jsd:jed)) ; sfc_state%internal_heat(:,:) = 0.0
364 if (alloc_iceshelves)
then
365 allocate(sfc_state%taux_shelf(isdb:iedb,jsd:jed)) ; sfc_state%taux_shelf(:,:) = 0.0
366 allocate(sfc_state%tauy_shelf(isd:ied,jsdb:jedb)) ; sfc_state%tauy_shelf(:,:) = 0.0
369 if (
present(gas_fields_ocn)) &
370 call coupler_type_spawn(gas_fields_ocn, sfc_state%tr_fields, &
371 (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
373 sfc_state%arrays_allocated = .true.
375 end subroutine allocate_surface_state
378 subroutine deallocate_surface_state(sfc_state)
381 if (.not.sfc_state%arrays_allocated)
return
383 if (
allocated(sfc_state%melt_potential))
deallocate(sfc_state%melt_potential)
384 if (
allocated(sfc_state%SST))
deallocate(sfc_state%SST)
385 if (
allocated(sfc_state%SSS))
deallocate(sfc_state%SSS)
386 if (
allocated(sfc_state%sfc_density))
deallocate(sfc_state%sfc_density)
387 if (
allocated(sfc_state%sea_lev))
deallocate(sfc_state%sea_lev)
388 if (
allocated(sfc_state%Hml))
deallocate(sfc_state%Hml)
389 if (
allocated(sfc_state%u))
deallocate(sfc_state%u)
390 if (
allocated(sfc_state%v))
deallocate(sfc_state%v)
391 if (
allocated(sfc_state%ocean_mass))
deallocate(sfc_state%ocean_mass)
392 if (
allocated(sfc_state%ocean_heat))
deallocate(sfc_state%ocean_heat)
393 if (
allocated(sfc_state%ocean_salt))
deallocate(sfc_state%ocean_salt)
394 if (
allocated(sfc_state%salt_deficit))
deallocate(sfc_state%salt_deficit)
396 call coupler_type_destructor(sfc_state%tr_fields)
398 sfc_state%arrays_allocated = .false.
400 end subroutine deallocate_surface_state
403 subroutine rotate_surface_state(sfc_state_in, G_in, sfc_state, G, turns)
406 type(
surface),
intent(inout) :: sfc_state
408 integer,
intent(in) :: turns
410 logical :: use_temperature, do_integrals, use_melt_potential, use_iceshelves
413 use_temperature =
allocated(sfc_state_in%SST) &
414 .and.
allocated(sfc_state_in%SSS)
415 use_melt_potential =
allocated(sfc_state_in%melt_potential)
416 do_integrals =
allocated(sfc_state_in%ocean_mass)
417 use_iceshelves =
allocated(sfc_state_in%taux_shelf) &
418 .and.
allocated(sfc_state_in%tauy_shelf)
420 if (.not. sfc_state%arrays_allocated)
then
421 call allocate_surface_state(sfc_state, g, &
422 use_temperature=use_temperature, &
423 do_integrals=do_integrals, &
424 use_meltpot=use_melt_potential, &
425 use_iceshelves=use_iceshelves &
427 sfc_state%arrays_allocated = .true.
430 if (use_temperature)
then
431 call rotate_array(sfc_state_in%SST, turns, sfc_state%SST)
432 call rotate_array(sfc_state_in%SSS, turns, sfc_state%SSS)
434 call rotate_array(sfc_state_in%sfc_density, turns, sfc_state%sfc_density)
437 call rotate_array(sfc_state_in%Hml, turns, sfc_state%Hml)
439 sfc_state%u, sfc_state%v)
440 call rotate_array(sfc_state_in%sea_lev, turns, sfc_state%sea_lev)
442 if (use_melt_potential)
then
443 call rotate_array(sfc_state_in%melt_potential, turns, sfc_state%melt_potential)
446 if (do_integrals)
then
447 call rotate_array(sfc_state_in%ocean_mass, turns, sfc_state%ocean_mass)
448 if (use_temperature)
then
449 call rotate_array(sfc_state_in%ocean_heat, turns, sfc_state%ocean_heat)
450 call rotate_array(sfc_state_in%ocean_salt, turns, sfc_state%ocean_salt)
451 call rotate_array(sfc_state_in%SSS, turns, sfc_state%TempxPmE)
452 call rotate_array(sfc_state_in%salt_deficit, turns, sfc_state%salt_deficit)
453 call rotate_array(sfc_state_in%internal_heat, turns, sfc_state%internal_heat)
457 if (use_iceshelves)
then
458 call rotate_vector(sfc_state_in%taux_shelf, sfc_state_in%tauy_shelf, turns, &
459 sfc_state%taux_shelf, sfc_state%tauy_shelf)
462 if (use_temperature .and.
allocated(sfc_state_in%frazil)) &
463 call rotate_array(sfc_state_in%frazil, turns, sfc_state%frazil)
466 sfc_state%T_is_conT = sfc_state_in%T_is_conT
467 sfc_state%S_is_absS = sfc_state_in%S_is_absS
470 if (coupler_type_initialized(sfc_state_in%tr_fields)) &
471 call mom_error(fatal,
"Rotation of surface state tracers is not yet " &
473 end subroutine rotate_surface_state
476 subroutine alloc_bt_cont_type(BT_cont, G, alloc_faces)
479 logical,
optional,
intent(in) :: alloc_faces
482 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
483 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
484 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
486 if (
associated(bt_cont))
call mom_error(fatal, &
487 "alloc_BT_cont_type called with an associated BT_cont_type pointer.")
490 allocate(bt_cont%FA_u_WW(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_WW(:,:) = 0.0
491 allocate(bt_cont%FA_u_W0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_W0(:,:) = 0.0
492 allocate(bt_cont%FA_u_E0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_E0(:,:) = 0.0
493 allocate(bt_cont%FA_u_EE(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_EE(:,:) = 0.0
494 allocate(bt_cont%uBT_WW(isdb:iedb,jsd:jed)) ; bt_cont%uBT_WW(:,:) = 0.0
495 allocate(bt_cont%uBT_EE(isdb:iedb,jsd:jed)) ; bt_cont%uBT_EE(:,:) = 0.0
497 allocate(bt_cont%FA_v_SS(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_SS(:,:) = 0.0
498 allocate(bt_cont%FA_v_S0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_S0(:,:) = 0.0
499 allocate(bt_cont%FA_v_N0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_N0(:,:) = 0.0
500 allocate(bt_cont%FA_v_NN(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_NN(:,:) = 0.0
501 allocate(bt_cont%vBT_SS(isd:ied,jsdb:jedb)) ; bt_cont%vBT_SS(:,:) = 0.0
502 allocate(bt_cont%vBT_NN(isd:ied,jsdb:jedb)) ; bt_cont%vBT_NN(:,:) = 0.0
504 if (
present(alloc_faces))
then ;
if (alloc_faces)
then
505 allocate(bt_cont%h_u(isdb:iedb,jsd:jed,1:g%ke)) ; bt_cont%h_u(:,:,:) = 0.0
506 allocate(bt_cont%h_v(isd:ied,jsdb:jedb,1:g%ke)) ; bt_cont%h_v(:,:,:) = 0.0
509 end subroutine alloc_bt_cont_type
512 subroutine dealloc_bt_cont_type(BT_cont)
515 if (.not.
associated(bt_cont))
return
517 deallocate(bt_cont%FA_u_WW) ;
deallocate(bt_cont%FA_u_W0)
518 deallocate(bt_cont%FA_u_E0) ;
deallocate(bt_cont%FA_u_EE)
519 deallocate(bt_cont%uBT_WW) ;
deallocate(bt_cont%uBT_EE)
521 deallocate(bt_cont%FA_v_SS) ;
deallocate(bt_cont%FA_v_S0)
522 deallocate(bt_cont%FA_v_N0) ;
deallocate(bt_cont%FA_v_NN)
523 deallocate(bt_cont%vBT_SS) ;
deallocate(bt_cont%vBT_NN)
525 if (
allocated(bt_cont%h_u))
deallocate(bt_cont%h_u)
526 if (
allocated(bt_cont%h_v))
deallocate(bt_cont%h_v)
530 end subroutine dealloc_bt_cont_type
533 subroutine mom_thermovar_chksum(mesg, tv, G)
534 character(len=*),
intent(in) :: mesg
541 if (
associated(tv%T)) &
542 call hchksum(tv%T, mesg//
" tv%T", g%HI)
543 if (
associated(tv%S)) &
544 call hchksum(tv%S, mesg//
" tv%S", g%HI)
545 if (
associated(tv%frazil)) &
546 call hchksum(tv%frazil, mesg//
" tv%frazil", g%HI, scale=g%US%Q_to_J_kg*g%US%RZ_to_kg_m2)
547 if (
associated(tv%salt_deficit)) &
548 call hchksum(tv%salt_deficit, mesg//
" tv%salt_deficit", g%HI, scale=g%US%RZ_to_kg_m2)
549 if (
associated(tv%TempxPmE)) &
550 call hchksum(tv%TempxPmE, mesg//
" tv%TempxPmE", g%HI, scale=g%US%RZ_to_kg_m2)
551 end subroutine mom_thermovar_chksum