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_extra_t => null()
247 real,
pointer,
dimension(:,:,:) :: kd_extra_s => null()
253 real,
pointer,
dimension(:,:,:) :: kd_shear => null()
256 real,
pointer,
dimension(:,:,:) :: kv_shear => null()
259 real,
pointer,
dimension(:,:,:) :: kv_shear_bu => null()
262 real,
pointer,
dimension(:,:,:) :: kv_slow => null()
265 real,
pointer,
dimension(:,:,:) :: tke_turb => null()
273 real,
allocatable :: fa_u_ee(:,:)
275 real,
allocatable :: fa_u_e0(:,:)
277 real,
allocatable :: fa_u_w0(:,:)
279 real,
allocatable :: fa_u_ww(:,:)
281 real,
allocatable :: ubt_ww(:,:)
283 real,
allocatable :: ubt_ee(:,:)
285 real,
allocatable :: fa_v_nn(:,:)
287 real,
allocatable :: fa_v_n0(:,:)
289 real,
allocatable :: fa_v_s0(:,:)
291 real,
allocatable :: fa_v_ss(:,:)
293 real,
allocatable :: vbt_ss(:,:)
295 real,
allocatable :: vbt_nn(:,:)
297 real,
allocatable :: h_u(:,:,:)
298 real,
allocatable :: h_v(:,:,:)
299 type(group_pass_type) :: pass_polarity_bt
300 type(group_pass_type) :: pass_fa_uv
307 subroutine allocate_surface_state(sfc_state, G, use_temperature, do_integrals, &
308 gas_fields_ocn, use_meltpot, use_iceshelves, omit_frazil)
310 type(
surface),
intent(inout) :: sfc_state
311 logical,
optional,
intent(in) :: use_temperature
312 logical,
optional,
intent(in) :: do_integrals
314 type(coupler_1d_bc_type), &
315 optional,
intent(in) :: gas_fields_ocn
320 logical,
optional,
intent(in) :: use_meltpot
321 logical,
optional,
intent(in) :: use_iceshelves
323 logical,
optional,
intent(in) :: omit_frazil
327 logical :: use_temp, alloc_integ, use_melt_potential, alloc_iceshelves, alloc_frazil
328 integer :: is, ie, js, je, isd, ied, jsd, jed
329 integer :: isdB, iedB, jsdB, jedB
331 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
332 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
333 isdb = g%isdB ; iedb = g%iedB; jsdb = g%jsdB ; jedb = g%jedB
335 use_temp = .true. ;
if (
present(use_temperature)) use_temp = use_temperature
336 alloc_integ = .true. ;
if (
present(do_integrals)) alloc_integ = do_integrals
337 use_melt_potential = .false. ;
if (
present(use_meltpot)) use_melt_potential = use_meltpot
338 alloc_iceshelves = .false. ;
if (
present(use_iceshelves)) alloc_iceshelves = use_iceshelves
339 alloc_frazil = .true. ;
if (
present(omit_frazil)) alloc_frazil = .not.omit_frazil
341 if (sfc_state%arrays_allocated)
return 344 allocate(sfc_state%SST(isd:ied,jsd:jed)) ; sfc_state%SST(:,:) = 0.0
345 allocate(sfc_state%SSS(isd:ied,jsd:jed)) ; sfc_state%SSS(:,:) = 0.0
347 allocate(sfc_state%sfc_density(isd:ied,jsd:jed)) ; sfc_state%sfc_density(:,:) = 0.0
349 if (use_temp .and. alloc_frazil)
then 350 allocate(sfc_state%frazil(isd:ied,jsd:jed)) ; sfc_state%frazil(:,:) = 0.0
352 allocate(sfc_state%sea_lev(isd:ied,jsd:jed)) ; sfc_state%sea_lev(:,:) = 0.0
353 allocate(sfc_state%Hml(isd:ied,jsd:jed)) ; sfc_state%Hml(:,:) = 0.0
354 allocate(sfc_state%u(isdb:iedb,jsd:jed)) ; sfc_state%u(:,:) = 0.0
355 allocate(sfc_state%v(isd:ied,jsdb:jedb)) ; sfc_state%v(:,:) = 0.0
357 if (use_melt_potential)
then 358 allocate(sfc_state%melt_potential(isd:ied,jsd:jed)) ; sfc_state%melt_potential(:,:) = 0.0
361 if (alloc_integ)
then 363 allocate(sfc_state%ocean_mass(isd:ied,jsd:jed)) ; sfc_state%ocean_mass(:,:) = 0.0
365 allocate(sfc_state%ocean_heat(isd:ied,jsd:jed)) ; sfc_state%ocean_heat(:,:) = 0.0
366 allocate(sfc_state%ocean_salt(isd:ied,jsd:jed)) ; sfc_state%ocean_salt(:,:) = 0.0
367 allocate(sfc_state%TempxPmE(isd:ied,jsd:jed)) ; sfc_state%TempxPmE(:,:) = 0.0
368 allocate(sfc_state%salt_deficit(isd:ied,jsd:jed)) ; sfc_state%salt_deficit(:,:) = 0.0
369 allocate(sfc_state%internal_heat(isd:ied,jsd:jed)) ; sfc_state%internal_heat(:,:) = 0.0
373 if (alloc_iceshelves)
then 374 allocate(sfc_state%taux_shelf(isdb:iedb,jsd:jed)) ; sfc_state%taux_shelf(:,:) = 0.0
375 allocate(sfc_state%tauy_shelf(isd:ied,jsdb:jedb)) ; sfc_state%tauy_shelf(:,:) = 0.0
378 if (
present(gas_fields_ocn)) &
379 call coupler_type_spawn(gas_fields_ocn, sfc_state%tr_fields, &
380 (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
382 sfc_state%arrays_allocated = .true.
384 end subroutine allocate_surface_state
387 subroutine deallocate_surface_state(sfc_state)
390 if (.not.sfc_state%arrays_allocated)
return 392 if (
allocated(sfc_state%melt_potential))
deallocate(sfc_state%melt_potential)
393 if (
allocated(sfc_state%SST))
deallocate(sfc_state%SST)
394 if (
allocated(sfc_state%SSS))
deallocate(sfc_state%SSS)
395 if (
allocated(sfc_state%sfc_density))
deallocate(sfc_state%sfc_density)
396 if (
allocated(sfc_state%sea_lev))
deallocate(sfc_state%sea_lev)
397 if (
allocated(sfc_state%Hml))
deallocate(sfc_state%Hml)
398 if (
allocated(sfc_state%u))
deallocate(sfc_state%u)
399 if (
allocated(sfc_state%v))
deallocate(sfc_state%v)
400 if (
allocated(sfc_state%ocean_mass))
deallocate(sfc_state%ocean_mass)
401 if (
allocated(sfc_state%ocean_heat))
deallocate(sfc_state%ocean_heat)
402 if (
allocated(sfc_state%ocean_salt))
deallocate(sfc_state%ocean_salt)
403 if (
allocated(sfc_state%salt_deficit))
deallocate(sfc_state%salt_deficit)
405 call coupler_type_destructor(sfc_state%tr_fields)
407 sfc_state%arrays_allocated = .false.
409 end subroutine deallocate_surface_state
412 subroutine rotate_surface_state(sfc_state_in, G_in, sfc_state, G, turns)
415 type(
surface),
intent(inout) :: sfc_state
417 integer,
intent(in) :: turns
419 logical :: use_temperature, do_integrals, use_melt_potential, use_iceshelves
422 use_temperature =
allocated(sfc_state_in%SST) &
423 .and.
allocated(sfc_state_in%SSS)
424 use_melt_potential =
allocated(sfc_state_in%melt_potential)
425 do_integrals =
allocated(sfc_state_in%ocean_mass)
426 use_iceshelves =
allocated(sfc_state_in%taux_shelf) &
427 .and.
allocated(sfc_state_in%tauy_shelf)
429 if (.not. sfc_state%arrays_allocated)
then 430 call allocate_surface_state(sfc_state, g, &
431 use_temperature=use_temperature, &
432 do_integrals=do_integrals, &
433 use_meltpot=use_melt_potential, &
434 use_iceshelves=use_iceshelves &
436 sfc_state%arrays_allocated = .true.
439 if (use_temperature)
then 440 call rotate_array(sfc_state_in%SST, turns, sfc_state%SST)
441 call rotate_array(sfc_state_in%SSS, turns, sfc_state%SSS)
443 call rotate_array(sfc_state_in%sfc_density, turns, sfc_state%sfc_density)
446 call rotate_array(sfc_state_in%Hml, turns, sfc_state%Hml)
448 sfc_state%u, sfc_state%v)
449 call rotate_array(sfc_state_in%sea_lev, turns, sfc_state%sea_lev)
451 if (use_melt_potential)
then 452 call rotate_array(sfc_state_in%melt_potential, turns, sfc_state%melt_potential)
455 if (do_integrals)
then 456 call rotate_array(sfc_state_in%ocean_mass, turns, sfc_state%ocean_mass)
457 if (use_temperature)
then 458 call rotate_array(sfc_state_in%ocean_heat, turns, sfc_state%ocean_heat)
459 call rotate_array(sfc_state_in%ocean_salt, turns, sfc_state%ocean_salt)
460 call rotate_array(sfc_state_in%SSS, turns, sfc_state%TempxPmE)
461 call rotate_array(sfc_state_in%salt_deficit, turns, sfc_state%salt_deficit)
462 call rotate_array(sfc_state_in%internal_heat, turns, sfc_state%internal_heat)
466 if (use_iceshelves)
then 467 call rotate_vector(sfc_state_in%taux_shelf, sfc_state_in%tauy_shelf, turns, &
468 sfc_state%taux_shelf, sfc_state%tauy_shelf)
471 if (use_temperature .and.
allocated(sfc_state_in%frazil)) &
472 call rotate_array(sfc_state_in%frazil, turns, sfc_state%frazil)
475 sfc_state%T_is_conT = sfc_state_in%T_is_conT
476 sfc_state%S_is_absS = sfc_state_in%S_is_absS
479 if (coupler_type_initialized(sfc_state_in%tr_fields)) &
480 call mom_error(fatal,
"Rotation of surface state tracers is not yet " &
482 end subroutine rotate_surface_state
485 subroutine alloc_bt_cont_type(BT_cont, G, alloc_faces)
488 logical,
optional,
intent(in) :: alloc_faces
491 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
492 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
493 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
495 if (
associated(bt_cont))
call mom_error(fatal, &
496 "alloc_BT_cont_type called with an associated BT_cont_type pointer.")
499 allocate(bt_cont%FA_u_WW(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_WW(:,:) = 0.0
500 allocate(bt_cont%FA_u_W0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_W0(:,:) = 0.0
501 allocate(bt_cont%FA_u_E0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_E0(:,:) = 0.0
502 allocate(bt_cont%FA_u_EE(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_EE(:,:) = 0.0
503 allocate(bt_cont%uBT_WW(isdb:iedb,jsd:jed)) ; bt_cont%uBT_WW(:,:) = 0.0
504 allocate(bt_cont%uBT_EE(isdb:iedb,jsd:jed)) ; bt_cont%uBT_EE(:,:) = 0.0
506 allocate(bt_cont%FA_v_SS(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_SS(:,:) = 0.0
507 allocate(bt_cont%FA_v_S0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_S0(:,:) = 0.0
508 allocate(bt_cont%FA_v_N0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_N0(:,:) = 0.0
509 allocate(bt_cont%FA_v_NN(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_NN(:,:) = 0.0
510 allocate(bt_cont%vBT_SS(isd:ied,jsdb:jedb)) ; bt_cont%vBT_SS(:,:) = 0.0
511 allocate(bt_cont%vBT_NN(isd:ied,jsdb:jedb)) ; bt_cont%vBT_NN(:,:) = 0.0
513 if (
present(alloc_faces))
then ;
if (alloc_faces)
then 514 allocate(bt_cont%h_u(isdb:iedb,jsd:jed,1:g%ke)) ; bt_cont%h_u(:,:,:) = 0.0
515 allocate(bt_cont%h_v(isd:ied,jsdb:jedb,1:g%ke)) ; bt_cont%h_v(:,:,:) = 0.0
518 end subroutine alloc_bt_cont_type
521 subroutine dealloc_bt_cont_type(BT_cont)
524 if (.not.
associated(bt_cont))
return 526 deallocate(bt_cont%FA_u_WW) ;
deallocate(bt_cont%FA_u_W0)
527 deallocate(bt_cont%FA_u_E0) ;
deallocate(bt_cont%FA_u_EE)
528 deallocate(bt_cont%uBT_WW) ;
deallocate(bt_cont%uBT_EE)
530 deallocate(bt_cont%FA_v_SS) ;
deallocate(bt_cont%FA_v_S0)
531 deallocate(bt_cont%FA_v_N0) ;
deallocate(bt_cont%FA_v_NN)
532 deallocate(bt_cont%vBT_SS) ;
deallocate(bt_cont%vBT_NN)
534 if (
allocated(bt_cont%h_u))
deallocate(bt_cont%h_u)
535 if (
allocated(bt_cont%h_v))
deallocate(bt_cont%h_v)
539 end subroutine dealloc_bt_cont_type
542 subroutine mom_thermovar_chksum(mesg, tv, G)
543 character(len=*),
intent(in) :: mesg
550 if (
associated(tv%T)) &
551 call hchksum(tv%T, mesg//
" tv%T", g%HI)
552 if (
associated(tv%S)) &
553 call hchksum(tv%S, mesg//
" tv%S", g%HI)
554 if (
associated(tv%frazil)) &
555 call hchksum(tv%frazil, mesg//
" tv%frazil", g%HI, scale=g%US%Q_to_J_kg*g%US%RZ_to_kg_m2)
556 if (
associated(tv%salt_deficit)) &
557 call hchksum(tv%salt_deficit, mesg//
" tv%salt_deficit", g%HI, scale=g%US%RZ_to_kg_m2)
558 if (
associated(tv%TempxPmE)) &
559 call hchksum(tv%TempxPmE, mesg//
" tv%TempxPmE", g%HI, scale=g%US%RZ_to_kg_m2)
560 end subroutine mom_thermovar_chksum
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
Vertical viscosities, drag coefficients, and related fields.
Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances.
A structure for creating arrays of pointers to 2D arrays.
Container for information about the summed layer transports and how they will vary as the barotropic ...
Describes the decomposed MOM domain and has routines for communications across PEs.
Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances.
Routines for error handling and I/O management.
Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90.
Provides subroutines for quantities specific to the equation of state.
The MOM_domain_type contains information about the domain decompositoin.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Provides transparent structures with groups of MOM6 variables and supporting routines.
Provides checksumming functions for debugging.
A control structure for the equation of state.
A structure for creating arrays of pointers to 3D arrays.