MOM6
ocean_model_mod Module Reference

Detailed Description

Top-level module for the MOM6 ocean model in coupled mode.

Data Types

interface  ocean_model_data_get
 This interface extracts a named scalar field or array from the ocean surface or public type. More...
 
type  ocean_public_type
 This type is used for communication with other components via the FMS coupler. The element names and types can be changed only with great deliberation, hence the persistnce of things like the cutsy element name "avg_kount". More...
 
type  ocean_state_type
 The ocean_state_type contains all information about the state of the ocean, with a format that is private so it can be readily changed without disrupting other coupled components. More...
 

Functions/Subroutines

subroutine, public ocean_model_init (Ocean_sfc, OS, Time_init, Time_in, wind_stagger, gas_fields_ocn)
 ocean_model_init initializes the ocean model, including registering fields for restarts and reading restart files if appropriate. More...
 
subroutine, public update_ocean_model (Ice_ocean_boundary, OS, Ocean_sfc, time_start_update, Ocean_coupling_time_step, update_dyn, update_thermo, Ocn_fluxes_used, start_cycle, end_cycle, cycle_length)
 update_ocean_model uses the forcing in Ice_ocean_boundary to advance the ocean model's state from the input value of Ocean_state (which must be for time time_start_update) for a time interval of Ocean_coupling_time_step, returning the publicly visible ocean surface properties in Ocean_sfc and storing the new ocean properties in Ocean_state. More...
 
subroutine, public ocean_model_restart (OS, timestamp)
 This subroutine writes out the ocean model restart file. More...
 
subroutine, public ocean_model_end (Ocean_sfc, Ocean_state, Time)
 ocean_model_end terminates the model run, saving the ocean state in a restart and deallocating any data associated with the ocean. More...
 
subroutine, public ocean_model_save_restart (OS, Time, directory, filename_suffix)
 ocean_model_save_restart causes restart files associated with the ocean to be written out. More...
 
subroutine initialize_ocean_public_type (input_domain, Ocean_sfc, diag, maskmap, gas_fields_ocn)
 Initialize the public ocean type. More...
 
subroutine convert_state_to_ocean_type (sfc_state, Ocean_sfc, G, US, patm, press_to_z)
 This subroutine translates the coupler's ocean_data_type into MOM's surface state variable. This may eventually be folded into the MOM code that calculates the surface state in the first place. Note the offset in the arrays because the ocean_data_type has no halo points in its arrays and always uses absolute indicies. More...
 
subroutine, public ocean_model_init_sfc (OS, Ocean_sfc)
 This subroutine extracts the surface properties from the ocean's internal state and stores them in the ocean type returned to the calling ice model. It has to be separate from the ocean_initialization call because the coupler module allocates the space for some of these variables. More...
 
subroutine, public ocean_model_flux_init (OS, verbosity)
 ocean_model_flux_init is used to initialize properties of the air-sea fluxes as determined by various run-time parameters. It can be called from non-ocean PEs, or PEs that have not yet been initialzed, and it can safely be called multiple times. More...
 
subroutine, public ocean_stock_pe (OS, index, value, time_index)
 Ocean_stock_pe - returns the integrated stocks of heat, water, etc. for conservation checks. Because of the way FMS is coded, only the root PE has the integrated amount, while all other PEs get 0. More...
 
subroutine ocean_model_data2d_get (OS, Ocean, name, array2D, isc, jsc)
 This subroutine extracts a named 2-D field from the ocean surface or public type. More...
 
subroutine ocean_model_data1d_get (OS, Ocean, name, value)
 This subroutine extracts a named scalar field from the ocean surface or public type. More...
 
subroutine, public ocean_public_type_chksum (id, timestep, ocn)
 Write out FMS-format checsums on fields from the ocean surface state. More...
 
subroutine, public get_ocean_grid (OS, Gridp)
 This subroutine gives a handle to the grid from ocean state. More...
 
subroutine, public ocean_model_get_uv_surf (OS, Ocean, name, array2D, isc, jsc)
 This subroutine extracts a named (u- or v-) 2-D surface current from ocean internal state. More...
 

Function/Subroutine Documentation

◆ convert_state_to_ocean_type()

subroutine ocean_model_mod::convert_state_to_ocean_type ( type(surface), intent(inout)  sfc_state,
type(ocean_public_type), intent(inout), target  Ocean_sfc,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
real, dimension(:,:), intent(in), optional  patm,
real, intent(in), optional  press_to_z 
)
private

This subroutine translates the coupler's ocean_data_type into MOM's surface state variable. This may eventually be folded into the MOM code that calculates the surface state in the first place. Note the offset in the arrays because the ocean_data_type has no halo points in its arrays and always uses absolute indicies.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]ocean_sfcA structure containing various publicly
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
[in]patmThe pressure at the ocean surface [Pa].
[in]press_to_zA conversion factor between pressure and ocean depth in m, usually 1/(rho_0*g) [m Pa-1].

Definition at line 821 of file ocean_model_MOM.F90.

821  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
822  !! describe the surface state of the ocean.
823  type(ocean_public_type), &
824  target, intent(inout) :: Ocean_sfc !< A structure containing various publicly
825  !! visible ocean surface fields, whose elements
826  !! have their data set here.
827  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
828  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
829  real, optional, intent(in) :: patm(:,:) !< The pressure at the ocean surface [Pa].
830  real, optional, intent(in) :: press_to_z !< A conversion factor between pressure and
831  !! ocean depth in m, usually 1/(rho_0*g) [m Pa-1].
832  ! Local variables
833  real :: IgR0
834  character(len=48) :: val_str
835  integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
836  integer :: i, j, i0, j0, is, ie, js, je
837 
838  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
839  call pass_vector(sfc_state%u, sfc_state%v, g%Domain)
840 
841  call mpp_get_compute_domain(ocean_sfc%Domain, isc_bnd, iec_bnd, &
842  jsc_bnd, jec_bnd)
843  if (present(patm)) then
844  ! Check that the inidicies in patm are (isc_bnd:iec_bnd,jsc_bnd:jec_bnd).
845  if (.not.present(press_to_z)) call mom_error(fatal, &
846  'convert_state_to_ocean_type: press_to_z must be present if patm is.')
847  endif
848 
849  i0 = is - isc_bnd ; j0 = js - jsc_bnd
850  if (sfc_state%T_is_conT) then
851  ! Convert the surface T from conservative T to potential T.
852  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
853  ocean_sfc%t_surf(i,j) = gsw_pt_from_ct(sfc_state%SSS(i+i0,j+j0), &
854  sfc_state%SST(i+i0,j+j0)) + celsius_kelvin_offset
855  enddo ; enddo
856  else
857  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
858  ocean_sfc%t_surf(i,j) = sfc_state%SST(i+i0,j+j0) + celsius_kelvin_offset
859  enddo ; enddo
860  endif
861  if (sfc_state%S_is_absS) then
862  ! Convert the surface S from absolute salinity to practical salinity.
863  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
864  ocean_sfc%s_surf(i,j) = gsw_sp_from_sr(sfc_state%SSS(i+i0,j+j0))
865  enddo ; enddo
866  else
867  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
868  ocean_sfc%s_surf(i,j) = sfc_state%SSS(i+i0,j+j0)
869  enddo ; enddo
870  endif
871 
872  if (present(patm)) then
873  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
874  ocean_sfc%sea_lev(i,j) = us%Z_to_m * sfc_state%sea_lev(i+i0,j+j0) + patm(i,j) * press_to_z
875  ocean_sfc%area(i,j) = us%L_to_m**2 * g%areaT(i+i0,j+j0)
876  enddo ; enddo
877  else
878  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
879  ocean_sfc%sea_lev(i,j) = us%Z_to_m * sfc_state%sea_lev(i+i0,j+j0)
880  ocean_sfc%area(i,j) = us%L_to_m**2 * g%areaT(i+i0,j+j0)
881  enddo ; enddo
882  endif
883 
884  if (allocated(sfc_state%frazil)) then
885  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
886  ocean_sfc%frazil(i,j) = us%Q_to_J_kg*us%RZ_to_kg_m2 * sfc_state%frazil(i+i0,j+j0)
887  enddo ; enddo
888  endif
889 
890  if (ocean_sfc%stagger == agrid) then
891  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
892  ocean_sfc%u_surf(i,j) = g%mask2dT(i+i0,j+j0) * us%L_T_to_m_s * &
893  0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i-1+i0,j+j0))
894  ocean_sfc%v_surf(i,j) = g%mask2dT(i+i0,j+j0) * us%L_T_to_m_s * &
895  0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0,j-1+j0))
896  enddo ; enddo
897  elseif (ocean_sfc%stagger == bgrid_ne) then
898  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
899  ocean_sfc%u_surf(i,j) = g%mask2dBu(i+i0,j+j0) * us%L_T_to_m_s * &
900  0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i+i0,j+j0+1))
901  ocean_sfc%v_surf(i,j) = g%mask2dBu(i+i0,j+j0) * us%L_T_to_m_s * &
902  0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0+1,j+j0))
903  enddo ; enddo
904  elseif (ocean_sfc%stagger == cgrid_ne) then
905  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
906  ocean_sfc%u_surf(i,j) = g%mask2dCu(i+i0,j+j0)*us%L_T_to_m_s * sfc_state%u(i+i0,j+j0)
907  ocean_sfc%v_surf(i,j) = g%mask2dCv(i+i0,j+j0)*us%L_T_to_m_s * sfc_state%v(i+i0,j+j0)
908  enddo ; enddo
909  else
910  write(val_str, '(I8)') ocean_sfc%stagger
911  call mom_error(fatal, "convert_state_to_ocean_type: "//&
912  "Ocean_sfc%stagger has the unrecognized value of "//trim(val_str))
913  endif
914 
915  if (coupler_type_initialized(sfc_state%tr_fields)) then
916  if (.not.coupler_type_initialized(ocean_sfc%fields)) then
917  call mom_error(fatal, "convert_state_to_ocean_type: "//&
918  "Ocean_sfc%fields has not been initialized.")
919  endif
920  call coupler_type_copy_data(sfc_state%tr_fields, ocean_sfc%fields)
921  endif
922 

◆ get_ocean_grid()

subroutine, public ocean_model_mod::get_ocean_grid ( type(ocean_state_type OS,
type(ocean_grid_type), pointer  Gridp 
)

This subroutine gives a handle to the grid from ocean state.

Parameters
osA structure containing the internal ocean state
gridpThe ocean's grid structure

Definition at line 1120 of file ocean_model_MOM.F90.

1120  ! Obtain the ocean grid.
1121  type(ocean_state_type) :: OS !< A structure containing the
1122  !! internal ocean state
1123  type(ocean_grid_type) , pointer :: Gridp !< The ocean's grid structure
1124 
1125  gridp => os%grid
1126  return

◆ initialize_ocean_public_type()

subroutine ocean_model_mod::initialize_ocean_public_type ( type(domain2d), intent(in)  input_domain,
type(ocean_public_type), intent(inout)  Ocean_sfc,
type(diag_ctrl), intent(in)  diag,
logical, dimension(:,:), intent(in), optional  maskmap,
type(coupler_1d_bc_type), intent(in), optional  gas_fields_ocn 
)
private

Initialize the public ocean type.

Parameters
[in]input_domainThe ocean model domain description
[in,out]ocean_sfcA structure containing various publicly visible ocean surface properties after initialization, whose elements are allocated here.
[in]diagA structure that regulates diagnsotic output
[in]maskmapA mask indicating which virtual processors
[in]gas_fields_ocnIf present, this type describes the

Definition at line 763 of file ocean_model_MOM.F90.

763  type(domain2D), intent(in) :: input_domain !< The ocean model domain description
764  type(ocean_public_type), intent(inout) :: Ocean_sfc !< A structure containing various publicly
765  !! visible ocean surface properties after initialization, whose
766  !! elements are allocated here.
767  type(diag_ctrl), intent(in) :: diag !< A structure that regulates diagnsotic output
768  logical, dimension(:,:), &
769  optional, intent(in) :: maskmap !< A mask indicating which virtual processors
770  !! are actually in use. If missing, all are used.
771  type(coupler_1d_bc_type), &
772  optional, intent(in) :: gas_fields_ocn !< If present, this type describes the
773  !! ocean and surface-ice fields that will participate
774  !! in the calculation of additional gas or other
775  !! tracer fluxes.
776 
777  integer :: xsz, ysz, layout(2)
778  ! ice-ocean-boundary fields are always allocated using absolute indicies
779  ! and have no halos.
780  integer :: isc, iec, jsc, jec
781 
782  call mpp_get_layout(input_domain,layout)
783  call mpp_get_global_domain(input_domain, xsize=xsz, ysize=ysz)
784  if (PRESENT(maskmap)) then
785  call mpp_define_domains((/1,xsz,1,ysz/),layout,ocean_sfc%Domain, maskmap=maskmap)
786  else
787  call mpp_define_domains((/1,xsz,1,ysz/),layout,ocean_sfc%Domain)
788  endif
789  call mpp_get_compute_domain(ocean_sfc%Domain, isc, iec, jsc, jec)
790 
791  allocate ( ocean_sfc%t_surf (isc:iec,jsc:jec), &
792  ocean_sfc%s_surf (isc:iec,jsc:jec), &
793  ocean_sfc%u_surf (isc:iec,jsc:jec), &
794  ocean_sfc%v_surf (isc:iec,jsc:jec), &
795  ocean_sfc%sea_lev(isc:iec,jsc:jec), &
796  ocean_sfc%area (isc:iec,jsc:jec), &
797  ocean_sfc%frazil (isc:iec,jsc:jec))
798 
799  ocean_sfc%t_surf(:,:) = 0.0 ! time averaged sst (Kelvin) passed to atmosphere/ice model
800  ocean_sfc%s_surf(:,:) = 0.0 ! time averaged sss (psu) passed to atmosphere/ice models
801  ocean_sfc%u_surf(:,:) = 0.0 ! time averaged u-current (m/sec) passed to atmosphere/ice models
802  ocean_sfc%v_surf(:,:) = 0.0 ! time averaged v-current (m/sec) passed to atmosphere/ice models
803  ocean_sfc%sea_lev(:,:) = 0.0 ! time averaged thickness of top model grid cell (m) plus patm/rho0/grav
804  ocean_sfc%frazil(:,:) = 0.0 ! time accumulated frazil (J/m^2) passed to ice model
805  ocean_sfc%area(:,:) = 0.0
806  ocean_sfc%axes = diag%axesT1%handles !diag axes to be used by coupler tracer flux diagnostics
807 
808  if (present(gas_fields_ocn)) then
809  call coupler_type_spawn(gas_fields_ocn, ocean_sfc%fields, (/isc,isc,iec,iec/), &
810  (/jsc,jsc,jec,jec/), suffix = '_ocn', as_needed=.true.)
811  endif
812 

◆ ocean_model_data1d_get()

subroutine ocean_model_mod::ocean_model_data1d_get ( type(ocean_state_type), pointer  OS,
type(ocean_public_type), intent(in)  Ocean,
character(len=*), intent(in)  name,
real, intent(out)  value 
)

This subroutine extracts a named scalar field from the ocean surface or public type.

Parameters
osA pointer to the structure containing the internal ocean state (intent in).
[in]oceanA structure containing various publicly visible ocean surface fields.
[in]nameThe name of the field to extract
[out]valueThe value of the named field

Definition at line 1075 of file ocean_model_MOM.F90.

1075  type(ocean_state_type), pointer :: OS !< A pointer to the structure containing the
1076  !! internal ocean state (intent in).
1077  type(ocean_public_type), intent(in) :: Ocean !< A structure containing various publicly
1078  !! visible ocean surface fields.
1079  character(len=*) , intent(in) :: name !< The name of the field to extract
1080  real , intent(out):: value !< The value of the named field
1081 
1082  if (.not.associated(os)) return
1083  if (.not.os%is_ocean_pe) return
1084 
1085  select case(name)
1086  case('c_p')
1087  value = os%C_p
1088  case default
1089  call mom_error(fatal,'get_ocean_grid_data1D: unknown argument name='//name)
1090  end select
1091 

◆ ocean_model_data2d_get()

subroutine ocean_model_mod::ocean_model_data2d_get ( type(ocean_state_type), pointer  OS,
type(ocean_public_type), intent(in)  Ocean,
character(len=*), intent(in)  name,
real, dimension(isc:,jsc:), intent(out)  array2D,
integer, intent(in)  isc,
integer, intent(in)  jsc 
)

This subroutine extracts a named 2-D field from the ocean surface or public type.

Parameters
osA pointer to the structure containing the internal ocean state (intent in).
[in]oceanA structure containing various publicly visible ocean surface fields.
[in]nameThe name of the field to extract
[out]array2dThe values of the named field, it must cover only the computational domain
[in]iscThe starting i-index of array2D
[in]jscThe starting j-index of array2D

Definition at line 1013 of file ocean_model_MOM.F90.

1013  use mom_constants, only : celsius_kelvin_offset
1014  type(ocean_state_type), pointer :: OS !< A pointer to the structure containing the
1015  !! internal ocean state (intent in).
1016  type(ocean_public_type), intent(in) :: Ocean !< A structure containing various publicly
1017  !! visible ocean surface fields.
1018  character(len=*) , intent(in) :: name !< The name of the field to extract
1019  real, dimension(isc:,jsc:), intent(out):: array2D !< The values of the named field, it must
1020  !! cover only the computational domain
1021  integer , intent(in) :: isc !< The starting i-index of array2D
1022  integer , intent(in) :: jsc !< The starting j-index of array2D
1023 
1024  integer :: g_isc, g_iec, g_jsc, g_jec,g_isd, g_ied, g_jsd, g_jed, i, j
1025 
1026  if (.not.associated(os)) return
1027  if (.not.os%is_ocean_pe) return
1028 
1029 ! The problem is %areaT is on MOM domain but Ice_Ocean_Boundary%... is on mpp domain.
1030 ! We want to return the MOM data on the mpp (compute) domain
1031 ! Get MOM domain extents
1032  call mpp_get_compute_domain(os%grid%Domain%mpp_domain, g_isc, g_iec, g_jsc, g_jec)
1033  call mpp_get_data_domain (os%grid%Domain%mpp_domain, g_isd, g_ied, g_jsd, g_jed)
1034 
1035  g_isc = g_isc-g_isd+1 ; g_iec = g_iec-g_isd+1 ; g_jsc = g_jsc-g_jsd+1 ; g_jec = g_jec-g_jsd+1
1036 
1037 
1038  select case(name)
1039  case('area')
1040  array2d(isc:,jsc:) = os%US%L_to_m**2*os%grid%areaT(g_isc:g_iec,g_jsc:g_jec)
1041  case('mask')
1042  array2d(isc:,jsc:) = os%grid%mask2dT(g_isc:g_iec,g_jsc:g_jec)
1043 !OR same result
1044 ! do j=g_jsc,g_jec ; do i=g_isc,g_iec
1045 ! array2D(isc+i-g_isc,jsc+j-g_jsc) = OS%grid%mask2dT(i,j)
1046 ! enddo ; enddo
1047  case('t_surf')
1048  array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1049  case('t_pme')
1050  array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1051  case('t_runoff')
1052  array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1053  case('t_calving')
1054  array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1055  case('btfHeat')
1056  array2d(isc:,jsc:) = 0
1057  case('cos_rot')
1058  array2d(isc:,jsc:) = os%grid%cos_rot(g_isc:g_iec,g_jsc:g_jec) ! =1
1059  case('sin_rot')
1060  array2d(isc:,jsc:) = os%grid%sin_rot(g_isc:g_iec,g_jsc:g_jec) ! =0
1061  case('s_surf')
1062  array2d(isc:,jsc:) = ocean%s_surf(isc:,jsc:)
1063  case('sea_lev')
1064  array2d(isc:,jsc:) = ocean%sea_lev(isc:,jsc:)
1065  case('frazil')
1066  array2d(isc:,jsc:) = ocean%frazil(isc:,jsc:)
1067  case default
1068  call mom_error(fatal,'get_ocean_grid_data2D: unknown argument name='//name)
1069  end select
1070 
Provides a few physical constants.

◆ ocean_model_end()

subroutine, public ocean_model_mod::ocean_model_end ( type(ocean_public_type), intent(inout)  Ocean_sfc,
type(ocean_state_type), pointer  Ocean_state,
type(time_type), intent(in)  Time 
)

ocean_model_end terminates the model run, saving the ocean state in a restart and deallocating any data associated with the ocean.

Parameters
[in,out]ocean_sfcAn ocean_public_type structure that is to be deallocated upon termination.
ocean_stateA pointer to the structure containing the internal ocean state to be deallocated upon termination.
[in]timeThe model time, used for writing restarts.

Definition at line 711 of file ocean_model_MOM.F90.

711  type(ocean_public_type), intent(inout) :: Ocean_sfc !< An ocean_public_type structure that is
712  !! to be deallocated upon termination.
713  type(ocean_state_type), pointer :: Ocean_state !< A pointer to the structure containing
714  !! the internal ocean state to be deallocated
715  !! upon termination.
716  type(time_type), intent(in) :: Time !< The model time, used for writing restarts.
717 
718  call ocean_model_save_restart(ocean_state, time)
719  call diag_mediator_end(time, ocean_state%diag)
720  call mom_end(ocean_state%MOM_CSp)
721  if (ocean_state%use_ice_shelf) call ice_shelf_end(ocean_state%Ice_shelf_CSp)

◆ ocean_model_flux_init()

subroutine, public ocean_model_mod::ocean_model_flux_init ( type(ocean_state_type), optional, pointer  OS,
integer, intent(in), optional  verbosity 
)

ocean_model_flux_init is used to initialize properties of the air-sea fluxes as determined by various run-time parameters. It can be called from non-ocean PEs, or PEs that have not yet been initialzed, and it can safely be called multiple times.

Parameters
osAn optional pointer to the ocean state, used to figure out if this is an ocean PE that has already been initialized.
[in]verbosityA 0-9 integer indicating a level of verbosity.

Definition at line 951 of file ocean_model_MOM.F90.

951  type(ocean_state_type), optional, pointer :: OS !< An optional pointer to the ocean state,
952  !! used to figure out if this is an ocean PE that
953  !! has already been initialized.
954  integer, optional, intent(in) :: verbosity !< A 0-9 integer indicating a level of verbosity.
955 
956  logical :: OS_is_set
957  integer :: verbose
958 
959  os_is_set = .false. ; if (present(os)) os_is_set = associated(os)
960 
961  ! Use this to control the verbosity of output; consider rethinking this logic later.
962  verbose = 5 ; if (os_is_set) verbose = 3
963  if (present(verbosity)) verbose = verbosity
964 
965  call call_tracer_flux_init(verbosity=verbose)
966 

◆ ocean_model_get_uv_surf()

subroutine, public ocean_model_mod::ocean_model_get_uv_surf ( type(ocean_state_type), pointer  OS,
type(ocean_public_type), intent(in)  Ocean,
character(len=*), intent(in)  name,
real, dimension(isc:,jsc:), intent(out)  array2D,
integer, intent(in)  isc,
integer, intent(in)  jsc 
)

This subroutine extracts a named (u- or v-) 2-D surface current from ocean internal state.

Parameters
osA pointer to the structure containing the internal ocean state (intent in).
[in]oceanA structure containing various publicly visible ocean surface fields.
[in]nameThe name of the current (ua or va) to extract
[out]array2dThe values of the named field, it must cover only the computational domain
[in]iscThe starting i-index of array2D
[in]jscThe starting j-index of array2D

Definition at line 1131 of file ocean_model_MOM.F90.

1131 
1132  type(ocean_state_type), pointer :: OS !< A pointer to the structure containing the
1133  !! internal ocean state (intent in).
1134  type(ocean_public_type), intent(in) :: Ocean !< A structure containing various publicly
1135  !! visible ocean surface fields.
1136  character(len=*) , intent(in) :: name !< The name of the current (ua or va) to extract
1137  real, dimension(isc:,jsc:), intent(out):: array2D !< The values of the named field, it must
1138  !! cover only the computational domain
1139  integer , intent(in) :: isc !< The starting i-index of array2D
1140  integer , intent(in) :: jsc !< The starting j-index of array2D
1141 
1142  type(ocean_grid_type) , pointer :: G !< The ocean's grid structure
1143  type(surface), pointer :: sfc_state !< A structure containing fields that
1144  !! describe the surface state of the ocean.
1145 
1146  integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
1147  integer :: i, j, i0, j0
1148  integer :: is, ie, js, je
1149 
1150  if (.not.associated(os)) return
1151  if (.not.os%is_ocean_pe) return
1152 
1153  g => os%grid
1154  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1155 
1156  call mpp_get_compute_domain(ocean%Domain, isc_bnd, iec_bnd, &
1157  jsc_bnd, jec_bnd)
1158 
1159  i0 = is - isc_bnd ; j0 = js - jsc_bnd
1160 
1161  sfc_state => os%sfc_state
1162 
1163  select case(name)
1164  case('ua')
1165  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
1166  array2d(i,j) = g%mask2dT(i+i0,j+j0) * &
1167  0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i-1+i0,j+j0))
1168  enddo ; enddo
1169  case('va')
1170  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
1171  array2d(i,j) = g%mask2dT(i+i0,j+j0) * &
1172  0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0,j-1+j0))
1173  enddo ; enddo
1174  case('ub')
1175  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
1176  array2d(i,j) = g%mask2dBu(i+i0,j+j0) * &
1177  0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i+i0,j+j0+1))
1178  enddo ; enddo
1179  case('vb')
1180  do j=jsc_bnd,jec_bnd ; do i=isc_bnd,iec_bnd
1181  array2d(i,j) = g%mask2dBu(i+i0,j+j0) * &
1182  0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0+1,j+j0))
1183  enddo ; enddo
1184  case default
1185  call mom_error(fatal,'ocean_model_get_UV_surf: unknown argument name='//name)
1186  end select
1187 

◆ ocean_model_init()

subroutine, public ocean_model_mod::ocean_model_init ( type(ocean_public_type), intent(inout), target  Ocean_sfc,
type(ocean_state_type), pointer  OS,
type(time_type), intent(in)  Time_init,
type(time_type), intent(in)  Time_in,
integer, intent(in), optional  wind_stagger,
type(coupler_1d_bc_type), intent(in), optional  gas_fields_ocn 
)

ocean_model_init initializes the ocean model, including registering fields for restarts and reading restart files if appropriate.

This subroutine initializes both the ocean state and the ocean surface type. Because of the way that indicies and domains are handled, Ocean_sfc must have been used in a previous call to initialize_ocean_type.

Parameters
[in,out]ocean_sfcA structure containing various publicly
osA structure whose internal contents are private to ocean_model_mod that may be used to contain all information about the ocean's interior state.
[in]time_initThe start time for the coupled model's calendar
[in]time_inThe time at which to initialize the ocean model.
[in]wind_staggerIf present, the staggering of the winds that are being provided in calls to update_ocean_model
[in]gas_fields_ocnIf present, this type describes the

Definition at line 228 of file ocean_model_MOM.F90.

228  type(ocean_public_type), target, &
229  intent(inout) :: Ocean_sfc !< A structure containing various publicly
230  !! visible ocean surface properties after initialization,
231  !! the data in this type is intent out.
232  type(ocean_state_type), pointer :: OS !< A structure whose internal
233  !! contents are private to ocean_model_mod that may be used to
234  !! contain all information about the ocean's interior state.
235  type(time_type), intent(in) :: Time_init !< The start time for the coupled model's calendar
236  type(time_type), intent(in) :: Time_in !< The time at which to initialize the ocean model.
237  integer, optional, intent(in) :: wind_stagger !< If present, the staggering of the winds that are
238  !! being provided in calls to update_ocean_model
239  type(coupler_1d_bc_type), &
240  optional, intent(in) :: gas_fields_ocn !< If present, this type describes the
241  !! ocean and surface-ice fields that will participate
242  !! in the calculation of additional gas or other
243  !! tracer fluxes, and can be used to spawn related
244  !! internal variables in the ice model.
245  ! Local variables
246  real :: Rho0 ! The Boussinesq ocean density [kg m-3].
247  real :: G_Earth ! The gravitational acceleration [m s-2].
248  real :: HFrz !< If HFrz > 0 (m), melt potential will be computed.
249  !! The actual depth over which melt potential is computed will
250  !! min(HFrz, OBLD), where OBLD is the boundary layer depth.
251  !! If HFrz <= 0 (default), melt potential will not be computed.
252  logical :: use_melt_pot!< If true, allocate melt_potential array
253 
254 ! This include declares and sets the variable "version".
255 #include "version_variable.h"
256  character(len=40) :: mdl = "ocean_model_init" ! This module's name.
257  character(len=48) :: stagger ! A string indicating the staggering locations for the
258  ! surface velocities returned to the coupler.
259  type(param_file_type) :: param_file !< A structure to parse for run-time parameters
260  logical :: use_temperature ! If true, temperature and salinity are state variables.
261 
262  call calltree_enter("ocean_model_init(), ocean_model_MOM.F90")
263  if (associated(os)) then
264  call mom_error(warning, "ocean_model_init called with an associated "// &
265  "ocean_state_type structure. Model is already initialized.")
266  return
267  endif
268  allocate(os)
269 
270  os%is_ocean_pe = ocean_sfc%is_ocean_pe
271  if (.not.os%is_ocean_pe) return
272 
273  os%Time = time_in ; os%Time_dyn = time_in
274  call initialize_mom(os%Time, time_init, param_file, os%dirs, os%MOM_CSp, &
275  os%restart_CSp, time_in, offline_tracer_mode=os%offline_tracer_mode, &
276  diag_ptr=os%diag, count_calls=.true.)
277  call get_mom_state_elements(os%MOM_CSp, g=os%grid, gv=os%GV, us=os%US, c_p=os%C_p, &
278  c_p_scaled=os%fluxes%C_p, use_temp=use_temperature)
279 
280  ! Read all relevant parameters and write them to the model log.
281  call log_version(param_file, mdl, version, "")
282 
283  call get_param(param_file, mdl, "SINGLE_STEPPING_CALL", os%single_step_call, &
284  "If true, advance the state of MOM with a single step "//&
285  "including both dynamics and thermodynamics. If false, "//&
286  "the two phases are advanced with separate calls.", default=.true.)
287  call get_param(param_file, mdl, "DT", os%dt, &
288  "The (baroclinic) dynamics time step. The time-step that "//&
289  "is actually used will be an integer fraction of the "//&
290  "forcing time-step.", units="s", fail_if_missing=.true.)
291  call get_param(param_file, mdl, "DT_THERM", os%dt_therm, &
292  "The thermodynamic and tracer advection time step. "//&
293  "Ideally DT_THERM should be an integer multiple of DT "//&
294  "and less than the forcing or coupling time-step, unless "//&
295  "THERMO_SPANS_COUPLING is true, in which case DT_THERM "//&
296  "can be an integer multiple of the coupling timestep. By "//&
297  "default DT_THERM is set to DT.", units="s", default=os%dt)
298  call get_param(param_file, "MOM", "THERMO_SPANS_COUPLING", os%thermo_spans_coupling, &
299  "If true, the MOM will take thermodynamic and tracer "//&
300  "timesteps that can be longer than the coupling timestep. "//&
301  "The actual thermodynamic timestep that is used in this "//&
302  "case is the largest integer multiple of the coupling "//&
303  "timestep that is less than or equal to DT_THERM.", default=.false.)
304  call get_param(param_file, mdl, "DIABATIC_FIRST", os%diabatic_first, &
305  "If true, apply diabatic and thermodynamic processes, "//&
306  "including buoyancy forcing and mass gain or loss, "//&
307  "before stepping the dynamics forward.", default=.false.)
308 
309  call get_param(param_file, mdl, "RESTART_CONTROL", os%Restart_control, &
310  "An integer whose bits encode which restart files are "//&
311  "written. Add 2 (bit 1) for a time-stamped file, and odd "//&
312  "(bit 0) for a non-time-stamped file. A restart file "//&
313  "will be saved at the end of the run segment for any "//&
314  "non-negative value.", default=1)
315  call get_param(param_file, mdl, "OCEAN_SURFACE_STAGGER", stagger, &
316  "A case-insensitive character string to indicate the "//&
317  "staggering of the surface velocity field that is "//&
318  "returned to the coupler. Valid values include "//&
319  "'A', 'B', or 'C'.", default="C")
320  if (uppercase(stagger(1:1)) == 'A') then ; ocean_sfc%stagger = agrid
321  elseif (uppercase(stagger(1:1)) == 'B') then ; ocean_sfc%stagger = bgrid_ne
322  elseif (uppercase(stagger(1:1)) == 'C') then ; ocean_sfc%stagger = cgrid_ne
323  else ; call mom_error(fatal,"ocean_model_init: OCEAN_SURFACE_STAGGER = "// &
324  trim(stagger)//" is invalid.") ; endif
325 
326  call get_param(param_file, mdl, "RHO_0", rho0, &
327  "The mean ocean density used with BOUSSINESQ true to "//&
328  "calculate accelerations and the mass for conservation "//&
329  "properties, or with BOUSSINSEQ false to convert some "//&
330  "parameters from vertical units of m to kg m-2.", &
331  units="kg m-3", default=1035.0)
332  call get_param(param_file, mdl, "G_EARTH", g_earth, &
333  "The gravitational acceleration of the Earth.", &
334  units="m s-2", default = 9.80)
335 
336  call get_param(param_file, mdl, "ICE_SHELF", os%use_ice_shelf, &
337  "If true, enables the ice shelf model.", default=.false.)
338 
339  call get_param(param_file, mdl, "ICEBERGS_APPLY_RIGID_BOUNDARY", os%icebergs_alter_ocean, &
340  "If true, allows icebergs to change boundary condition felt by ocean", default=.false.)
341 
342  os%press_to_z = 1.0/(rho0*g_earth)
343 
344  ! Consider using a run-time flag to determine whether to do the diagnostic
345  ! vertical integrals, since the related 3-d sums are not negligible in cost.
346  call get_param(param_file, mdl, "HFREEZE", hfrz, &
347  "If HFREEZE > 0, melt potential will be computed. The actual depth "//&
348  "over which melt potential is computed will be min(HFREEZE, OBLD), "//&
349  "where OBLD is the boundary layer depth. If HFREEZE <= 0 (default), "//&
350  "melt potential will not be computed.", units="m", default=-1.0, do_not_log=.true.)
351 
352  if (hfrz .gt. 0.0) then
353  use_melt_pot=.true.
354  else
355  use_melt_pot=.false.
356  endif
357 
358  call allocate_surface_state(os%sfc_state, os%grid, use_temperature, do_integrals=.true., &
359  gas_fields_ocn=gas_fields_ocn, use_meltpot=use_melt_pot)
360 
361  if (present(wind_stagger)) then
362  call surface_forcing_init(time_in, os%grid, os%US, param_file, os%diag, &
363  os%forcing_CSp, wind_stagger)
364  else
365  call surface_forcing_init(time_in, os%grid, os%US, param_file, os%diag, &
366  os%forcing_CSp)
367  endif
368 
369  if (os%use_ice_shelf) then
370  call initialize_ice_shelf(param_file, os%grid, os%Time, os%ice_shelf_CSp, &
371  os%diag, os%forces, os%fluxes)
372  endif
373  if (os%icebergs_alter_ocean) then
374  call marine_ice_init(os%Time, os%grid, param_file, os%diag, os%marine_ice_CSp)
375  if (.not. os%use_ice_shelf) &
376  call allocate_forcing_type(os%grid, os%fluxes, shelf=.true.)
377  endif
378 
379  call get_param(param_file, mdl, "USE_WAVES", os%Use_Waves, &
380  "If true, enables surface wave modules.", default=.false.)
381  if (os%use_waves) then
382  call mom_wave_interface_init(os%Time, os%grid, os%GV, os%US, param_file, os%Waves, os%diag)
383  else
384  call mom_wave_interface_init_lite(param_file)
385  endif
386 
387  if (associated(os%grid%Domain%maskmap)) then
388  call initialize_ocean_public_type(os%grid%Domain%mpp_domain, ocean_sfc, &
389  os%diag, maskmap=os%grid%Domain%maskmap, &
390  gas_fields_ocn=gas_fields_ocn)
391  else
392  call initialize_ocean_public_type(os%grid%Domain%mpp_domain, ocean_sfc, &
393  os%diag, gas_fields_ocn=gas_fields_ocn)
394  endif
395 
396  ! This call can only occur here if the coupler_bc_type variables have been
397  ! initialized already using the information from gas_fields_ocn.
398  if (present(gas_fields_ocn)) then
399  call coupler_type_set_diags(ocean_sfc%fields, "ocean_sfc", &
400  ocean_sfc%axes(1:2), time_in)
401 
402  call extract_surface_state(os%MOM_CSp, os%sfc_state)
403 
404  call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
405 
406  endif
407 
408  call close_param_file(param_file)
409  call diag_mediator_close_registration(os%diag)
410 
411  call mom_mesg('==== Completed MOM6 Coupled Initialization ====', 2)
412 
413  call calltree_leave("ocean_model_init(")

◆ ocean_model_init_sfc()

subroutine, public ocean_model_mod::ocean_model_init_sfc ( type(ocean_state_type), pointer  OS,
type(ocean_public_type), intent(inout)  Ocean_sfc 
)

This subroutine extracts the surface properties from the ocean's internal state and stores them in the ocean type returned to the calling ice model. It has to be separate from the ocean_initialization call because the coupler module allocates the space for some of these variables.

Parameters
osThe structure with the complete ocean state
[in,out]ocean_sfcA structure containing various publicly visible ocean surface properties after initialization, whose elements have their data set here.

Definition at line 930 of file ocean_model_MOM.F90.

930  type(ocean_state_type), pointer :: OS !< The structure with the complete ocean state
931  type(ocean_public_type), intent(inout) :: Ocean_sfc !< A structure containing various publicly
932  !! visible ocean surface properties after initialization, whose
933  !! elements have their data set here.
934  integer :: is, ie, js, je
935 
936  is = os%grid%isc ; ie = os%grid%iec ; js = os%grid%jsc ; je = os%grid%jec
937  call coupler_type_spawn(ocean_sfc%fields, os%sfc_state%tr_fields, &
938  (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
939 
940  call extract_surface_state(os%MOM_CSp, os%sfc_state)
941 
942  call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
943 

◆ ocean_model_restart()

subroutine, public ocean_model_mod::ocean_model_restart ( type(ocean_state_type), pointer  OS,
character(len=*), intent(in), optional  timestamp 
)

This subroutine writes out the ocean model restart file.

Parameters
osA pointer to the structure containing the internal ocean state being saved to a restart file
[in]timestampAn optional timestamp string that should be prepended to the file name. (Currently this is unused.)

Definition at line 673 of file ocean_model_MOM.F90.

673  type(ocean_state_type), pointer :: OS !< A pointer to the structure containing the
674  !! internal ocean state being saved to a restart file
675  character(len=*), optional, intent(in) :: timestamp !< An optional timestamp string that should be
676  !! prepended to the file name. (Currently this is unused.)
677 
678  if (.not.mom_state_is_synchronized(os%MOM_CSp)) &
679  call mom_error(warning, "End of MOM_main reached with inconsistent "//&
680  "dynamics and advective times. Additional restart fields "//&
681  "that have not been coded yet would be required for reproducibility.")
682  if (.not.os%fluxes%fluxes_used) call mom_error(fatal, "ocean_model_restart "//&
683  "was called with unused buoyancy fluxes. For conservation, the ocean "//&
684  "restart files can only be created after the buoyancy forcing is applied.")
685 
686  if (btest(os%Restart_control,1)) then
687  call save_restart(os%dirs%restart_output_dir, os%Time, os%grid, &
688  os%restart_CSp, .true., gv=os%GV)
689  call forcing_save_restart(os%forcing_CSp, os%grid, os%Time, &
690  os%dirs%restart_output_dir, .true.)
691  if (os%use_ice_shelf) then
692  call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir, .true.)
693  endif
694  endif
695  if (btest(os%Restart_control,0)) then
696  call save_restart(os%dirs%restart_output_dir, os%Time, os%grid, &
697  os%restart_CSp, gv=os%GV)
698  call forcing_save_restart(os%forcing_CSp, os%grid, os%Time, &
699  os%dirs%restart_output_dir)
700  if (os%use_ice_shelf) then
701  call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir)
702  endif
703  endif
704 

◆ ocean_model_save_restart()

subroutine, public ocean_model_mod::ocean_model_save_restart ( type(ocean_state_type), pointer  OS,
type(time_type), intent(in)  Time,
character(len=*), intent(in), optional  directory,
character(len=*), intent(in), optional  filename_suffix 
)

ocean_model_save_restart causes restart files associated with the ocean to be written out.

Parameters
osA pointer to the structure containing the internal ocean state (in).
[in]timeThe model time at this call, needed for mpp_write calls.
[in]directoryAn optional directory into which to write these restart files.
[in]filename_suffixAn optional suffix (e.g., a time-stamp) to append to the restart file names.

Definition at line 727 of file ocean_model_MOM.F90.

727  type(ocean_state_type), pointer :: OS !< A pointer to the structure containing the
728  !! internal ocean state (in).
729  type(time_type), intent(in) :: Time !< The model time at this call, needed for mpp_write calls.
730  character(len=*), optional, intent(in) :: directory !< An optional directory into which to
731  !! write these restart files.
732  character(len=*), optional, intent(in) :: filename_suffix !< An optional suffix (e.g., a time-stamp)
733  !! to append to the restart file names.
734 ! Note: This is a new routine - it will need to exist for the new incremental
735 ! checkpointing. It will also be called by ocean_model_end, giving the same
736 ! restart behavior as now in FMS.
737  character(len=200) :: restart_dir
738 
739  if (.not.mom_state_is_synchronized(os%MOM_CSp)) &
740  call mom_error(warning, "ocean_model_save_restart called with inconsistent "//&
741  "dynamics and advective times. Additional restart fields "//&
742  "that have not been coded yet would be required for reproducibility.")
743  if (.not.os%fluxes%fluxes_used) call mom_error(fatal, "ocean_model_save_restart "//&
744  "was called with unused buoyancy fluxes. For conservation, the ocean "//&
745  "restart files can only be created after the buoyancy forcing is applied.")
746 
747  if (present(directory)) then ; restart_dir = directory
748  else ; restart_dir = os%dirs%restart_output_dir ; endif
749 
750  call save_restart(restart_dir, time, os%grid, os%restart_CSp, gv=os%GV)
751 
752  call forcing_save_restart(os%forcing_CSp, os%grid, time, restart_dir)
753 
754  if (os%use_ice_shelf) then
755  call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir)
756  endif
757 

◆ ocean_public_type_chksum()

subroutine, public ocean_model_mod::ocean_public_type_chksum ( character(len=*), intent(in)  id,
integer, intent(in)  timestep,
type(ocean_public_type), intent(in)  ocn 
)

Write out FMS-format checsums on fields from the ocean surface state.

Parameters
[in]idAn identifying string for this call
[in]timestepThe number of elapsed timesteps
[in]ocnA structure containing various publicly visible ocean surface fields.

Definition at line 1096 of file ocean_model_MOM.F90.

1096 
1097  character(len=*), intent(in) :: id !< An identifying string for this call
1098  integer, intent(in) :: timestep !< The number of elapsed timesteps
1099  type(ocean_public_type), intent(in) :: ocn !< A structure containing various publicly
1100  !! visible ocean surface fields.
1101  integer :: n, m, outunit
1102 
1103  outunit = stdout()
1104 
1105  write(outunit,*) "BEGIN CHECKSUM(ocean_type):: ", id, timestep
1106  write(outunit,100) 'ocean%t_surf ',mpp_chksum(ocn%t_surf )
1107  write(outunit,100) 'ocean%s_surf ',mpp_chksum(ocn%s_surf )
1108  write(outunit,100) 'ocean%u_surf ',mpp_chksum(ocn%u_surf )
1109  write(outunit,100) 'ocean%v_surf ',mpp_chksum(ocn%v_surf )
1110  write(outunit,100) 'ocean%sea_lev ',mpp_chksum(ocn%sea_lev)
1111  write(outunit,100) 'ocean%frazil ',mpp_chksum(ocn%frazil )
1112 
1113  call coupler_type_write_chksums(ocn%fields, outunit, 'ocean%')
1114 100 FORMAT(" CHECKSUM::",a20," = ",z20)
1115 

◆ ocean_stock_pe()

subroutine, public ocean_model_mod::ocean_stock_pe ( type(ocean_state_type), pointer  OS,
integer, intent(in)  index,
real, intent(out)  value,
integer, intent(in), optional  time_index 
)

Ocean_stock_pe - returns the integrated stocks of heat, water, etc. for conservation checks. Because of the way FMS is coded, only the root PE has the integrated amount, while all other PEs get 0.

Parameters
osA structure containing the internal ocean state. The data in OS is intent in.
[in]indexThe stock index for the quantity of interest.
[out]valueSum returned for the conservation quantity of interest.
[in]time_indexAn unused optional argument, present only for interfacial compatibility with other models.

Definition at line 973 of file ocean_model_MOM.F90.

973  use stock_constants_mod, only : istock_water, istock_heat,istock_salt
974  type(ocean_state_type), pointer :: OS !< A structure containing the internal ocean state.
975  !! The data in OS is intent in.
976  integer, intent(in) :: index !< The stock index for the quantity of interest.
977  real, intent(out) :: value !< Sum returned for the conservation quantity of interest.
978  integer, optional, intent(in) :: time_index !< An unused optional argument, present only for
979  !! interfacial compatibility with other models.
980 ! Arguments: OS - A structure containing the internal ocean state.
981 ! (in) index - Index of conservation quantity of interest.
982 ! (in) value - Sum returned for the conservation quantity of interest.
983 ! (in,opt) time_index - Index for time level to use if this is necessary.
984 
985  real :: salt
986 
987  value = 0.0
988  if (.not.associated(os)) return
989  if (.not.os%is_ocean_pe) return
990 
991  select case (index)
992  case (istock_water) ! Return the mass of fresh water in the ocean in kg.
993  if (os%GV%Boussinesq) then
994  call get_ocean_stocks(os%MOM_CSp, mass=value, on_pe_only=.true.)
995  else ! In non-Boussinesq mode, the mass of salt needs to be subtracted.
996  call get_ocean_stocks(os%MOM_CSp, mass=value, salt=salt, on_pe_only=.true.)
997  value = value - salt
998  endif
999  case (istock_heat) ! Return the heat content of the ocean in J.
1000  call get_ocean_stocks(os%MOM_CSp, heat=value, on_pe_only=.true.)
1001  case (istock_salt) ! Return the mass of the salt in the ocean in kg.
1002  call get_ocean_stocks(os%MOM_CSp, salt=value, on_pe_only=.true.)
1003  case default ; value = 0.0
1004  end select
1005  ! If the FMS coupler is changed so that Ocean_stock_PE is only called on
1006  ! ocean PEs, uncomment the following and eliminate the on_PE_only flags above.
1007  ! if (.not.is_root_pe()) value = 0.0
1008 

◆ update_ocean_model()

subroutine, public ocean_model_mod::update_ocean_model ( type(ice_ocean_boundary_type), intent(in)  Ice_ocean_boundary,
type(ocean_state_type), pointer  OS,
type(ocean_public_type), intent(inout)  Ocean_sfc,
type(time_type), intent(in)  time_start_update,
type(time_type), intent(in)  Ocean_coupling_time_step,
logical, intent(in), optional  update_dyn,
logical, intent(in), optional  update_thermo,
logical, intent(in), optional  Ocn_fluxes_used,
logical, intent(in), optional  start_cycle,
logical, intent(in), optional  end_cycle,
real, intent(in), optional  cycle_length 
)

update_ocean_model uses the forcing in Ice_ocean_boundary to advance the ocean model's state from the input value of Ocean_state (which must be for time time_start_update) for a time interval of Ocean_coupling_time_step, returning the publicly visible ocean surface properties in Ocean_sfc and storing the new ocean properties in Ocean_state.

Parameters
[in]ice_ocean_boundaryA structure containing the various
osA pointer to a private structure containing the
[in,out]ocean_sfcA structure containing all the publicly visible
[in]time_start_updateThe time at the beginning of the update step.
[in]ocean_coupling_time_stepThe amount of time over which to advance the ocean.
[in]update_dynIf present and false, do not do updates due to the ocean dynamics.
[in]update_thermoIf present and false, do not do updates due to the ocean thermodynamics or remapping.
[in]ocn_fluxes_usedIf present, this indicates whether the cumulative thermodynamic fluxes from the ocean, like frazil, have been used and should be reset.
[in]start_cycleThis indicates whether this call is to be treated as the first call to step_MOM in a time-stepping cycle; missing is like true.
[in]end_cycleThis indicates whether this call is to be treated as the last call to step_MOM in a time-stepping cycle; missing is like true.
[in]cycle_lengthThe duration of a coupled time stepping cycle [s].

Definition at line 424 of file ocean_model_MOM.F90.

424  type(ice_ocean_boundary_type), &
425  intent(in) :: Ice_ocean_boundary !< A structure containing the various
426  !! forcing fields coming from the ice and atmosphere.
427  type(ocean_state_type), &
428  pointer :: OS !< A pointer to a private structure containing the
429  !! internal ocean state.
430  type(ocean_public_type), &
431  intent(inout) :: Ocean_sfc !< A structure containing all the publicly visible
432  !! ocean surface fields after a coupling time step.
433  !! The data in this type is intent out.
434  type(time_type), intent(in) :: time_start_update !< The time at the beginning of the update step.
435  type(time_type), intent(in) :: Ocean_coupling_time_step !< The amount of time over which to
436  !! advance the ocean.
437  logical, optional, intent(in) :: update_dyn !< If present and false, do not do updates
438  !! due to the ocean dynamics.
439  logical, optional, intent(in) :: update_thermo !< If present and false, do not do updates
440  !! due to the ocean thermodynamics or remapping.
441  logical, optional, intent(in) :: Ocn_fluxes_used !< If present, this indicates whether the
442  !! cumulative thermodynamic fluxes from the ocean,
443  !! like frazil, have been used and should be reset.
444  logical, optional, intent(in) :: start_cycle !< This indicates whether this call is to be
445  !! treated as the first call to step_MOM in a
446  !! time-stepping cycle; missing is like true.
447  logical, optional, intent(in) :: end_cycle !< This indicates whether this call is to be
448  !! treated as the last call to step_MOM in a
449  !! time-stepping cycle; missing is like true.
450  real, optional, intent(in) :: cycle_length !< The duration of a coupled time stepping cycle [s].
451 
452  ! Local variables
453  type(time_type) :: Time_seg_start ! Stores the dynamic or thermodynamic ocean model time at the
454  ! start of this call to allow step_MOM to temporarily change the time
455  ! as seen by internal modules.
456  type(time_type) :: Time_thermo_start ! Stores the ocean model thermodynamics time at the start of
457  ! this call to allow step_MOM to temporarily change the time as seen by
458  ! internal modules.
459  type(time_type) :: Time1 ! The value of the ocean model's time at the start of a call to step_MOM.
460  integer :: index_bnds(4) ! The computational domain index bounds in the ice-ocean boundary type.
461  real :: weight ! Flux accumulation weight of the current fluxes.
462  real :: dt_coupling ! The coupling time step [s].
463  real :: dt_therm ! A limited and quantized version of OS%dt_therm [s].
464  real :: dt_dyn ! The dynamics time step [s].
465  real :: dtdia ! The diabatic time step [s].
466  real :: t_elapsed_seg ! The elapsed time in this update segment [s].
467  integer :: n ! The internal iteration counter.
468  integer :: nts ! The number of baroclinic dynamics time steps in a thermodynamic step.
469  integer :: n_max ! The number of calls to step_MOM dynamics in this call to update_ocean_model.
470  integer :: n_last_thermo ! The iteration number the last time thermodynamics were updated.
471  logical :: thermo_does_span_coupling ! If true, thermodynamic forcing spans multiple dynamic timesteps.
472  logical :: do_dyn ! If true, step the ocean dynamics and transport.
473  logical :: do_thermo ! If true, step the ocean thermodynamics.
474  logical :: step_thermo ! If true, take a thermodynamic step.
475  integer :: is, ie, js, je
476 
477  call calltree_enter("update_ocean_model(), ocean_model_MOM.F90")
478  dt_coupling = time_type_to_real(ocean_coupling_time_step)
479 
480  if (.not.associated(os)) then
481  call mom_error(fatal, "update_ocean_model called with an unassociated "// &
482  "ocean_state_type structure. ocean_model_init must be "// &
483  "called first to allocate this structure.")
484  return
485  endif
486 
487  do_dyn = .true. ; if (present(update_dyn)) do_dyn = update_dyn
488  do_thermo = .true. ; if (present(update_thermo)) do_thermo = update_thermo
489 
490  if (do_thermo .and. (time_start_update /= os%Time)) &
491  call mom_error(warning, "update_ocean_model: internal clock does not "//&
492  "agree with time_start_update argument.")
493  if (do_dyn .and. (time_start_update /= os%Time_dyn)) &
494  call mom_error(warning, "update_ocean_model: internal dynamics clock does not "//&
495  "agree with time_start_update argument.")
496 
497  if (.not.(do_dyn .or. do_thermo)) call mom_error(fatal, &
498  "update_ocean_model called without updating either dynamics or thermodynamics.")
499  if (do_dyn .and. do_thermo .and. (os%Time /= os%Time_dyn)) call mom_error(fatal, &
500  "update_ocean_model called to update both dynamics and thermodynamics with inconsistent clocks.")
501 
502  ! This is benign but not necessary if ocean_model_init_sfc was called or if
503  ! OS%sfc_state%tr_fields was spawned in ocean_model_init. Consider removing it.
504  is = os%grid%isc ; ie = os%grid%iec ; js = os%grid%jsc ; je = os%grid%jec
505  call coupler_type_spawn(ocean_sfc%fields, os%sfc_state%tr_fields, &
506  (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
507 
508  ! Translate Ice_ocean_boundary into fluxes and forces.
509  call mpp_get_compute_domain(ocean_sfc%Domain, index_bnds(1), index_bnds(2), &
510  index_bnds(3), index_bnds(4))
511 
512  if (do_dyn) then
513  call convert_iob_to_forces(ice_ocean_boundary, os%forces, index_bnds, os%Time_dyn, os%grid, os%US, &
514  os%forcing_CSp, dt_forcing=dt_coupling, reset_avg=os%fluxes%fluxes_used)
515  if (os%use_ice_shelf) &
516  call add_shelf_forces(os%grid, os%US, os%Ice_shelf_CSp, os%forces)
517  if (os%icebergs_alter_ocean) &
518  call iceberg_forces(os%grid, os%forces, os%use_ice_shelf, &
519  os%sfc_state, dt_coupling, os%marine_ice_CSp)
520  endif
521 
522  if (do_thermo) then
523  if (os%fluxes%fluxes_used) then
524  call convert_iob_to_fluxes(ice_ocean_boundary, os%fluxes, index_bnds, os%Time, dt_coupling, &
525  os%grid, os%US, os%forcing_CSp, os%sfc_state)
526 
527  ! Add ice shelf fluxes
528  if (os%use_ice_shelf) &
529  call shelf_calc_flux(os%sfc_state, os%fluxes, os%Time, dt_coupling, os%Ice_shelf_CSp)
530  if (os%icebergs_alter_ocean) &
531  call iceberg_fluxes(os%grid, os%US, os%fluxes, os%use_ice_shelf, &
532  os%sfc_state, dt_coupling, os%marine_ice_CSp)
533 
534 #ifdef _USE_GENERIC_TRACER
535  call enable_averaging(dt_coupling, os%Time + ocean_coupling_time_step, os%diag) !Is this needed?
536  call mom_generic_tracer_fluxes_accumulate(os%fluxes, 1.0) ! Here weight=1, so just store the current fluxes
537  call disable_averaging(os%diag)
538 #endif
539  else
540  ! The previous fluxes have not been used yet, so translate the input fluxes
541  ! into a temporary type and then accumulate them in about 20 lines.
542  os%flux_tmp%C_p = os%fluxes%C_p
543  call convert_iob_to_fluxes(ice_ocean_boundary, os%flux_tmp, index_bnds, os%Time, dt_coupling, &
544  os%grid, os%US, os%forcing_CSp, os%sfc_state)
545 
546  if (os%use_ice_shelf) &
547  call shelf_calc_flux(os%sfc_state, os%flux_tmp, os%Time, dt_coupling, os%Ice_shelf_CSp)
548  if (os%icebergs_alter_ocean) &
549  call iceberg_fluxes(os%grid, os%US, os%flux_tmp, os%use_ice_shelf, &
550  os%sfc_state, dt_coupling, os%marine_ice_CSp)
551 
552  call fluxes_accumulate(os%flux_tmp, os%fluxes, os%grid, weight)
553 #ifdef _USE_GENERIC_TRACER
554  ! Incorporate the current tracer fluxes into the running averages
555  call mom_generic_tracer_fluxes_accumulate(os%flux_tmp, weight)
556 #endif
557  endif
558  endif
559 
560  ! The net mass forcing is not currently used in the MOM6 dynamics solvers, so this is may be unnecessary.
561  if (do_dyn .and. associated(os%forces%net_mass_src) .and. .not.os%forces%net_mass_src_set) &
562  call get_net_mass_forcing(os%fluxes, os%grid, os%US, os%forces%net_mass_src)
563 
564  if (os%use_waves .and. do_thermo) then
565  ! For now, the waves are only updated on the thermodynamics steps, because that is where
566  ! the wave intensities are actually used to drive mixing. At some point, the wave updates
567  ! might also need to become a part of the ocean dynamics, according to B. Reichl.
568  call update_surface_waves(os%grid, os%GV, os%US, os%time, ocean_coupling_time_step, os%waves)
569  endif
570 
571  if ((os%nstep==0) .and. (os%nstep_thermo==0)) then ! This is the first call to update_ocean_model.
572  call finish_mom_initialization(os%Time, os%dirs, os%MOM_CSp, os%restart_CSp)
573  endif
574 
575  time_thermo_start = os%Time
576  time_seg_start = os%Time ; if (do_dyn) time_seg_start = os%Time_dyn
577  time1 = time_seg_start
578 
579  if (os%offline_tracer_mode .and. do_thermo) then
580  call step_offline(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp)
581  elseif ((.not.do_thermo) .or. (.not.do_dyn)) then
582  ! The call sequence is being orchestrated from outside of update_ocean_model.
583  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp, &
584  waves=os%Waves, do_dynamics=do_dyn, do_thermodynamics=do_thermo, &
585  start_cycle=start_cycle, end_cycle=end_cycle, cycle_length=cycle_length, &
586  reset_therm=ocn_fluxes_used)
587  elseif (os%single_step_call) then
588  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp, waves=os%Waves)
589  else ! Step both the dynamics and thermodynamics with separate calls.
590  n_max = 1 ; if (dt_coupling > os%dt) n_max = ceiling(dt_coupling/os%dt - 0.001)
591  dt_dyn = dt_coupling / real(n_max)
592  thermo_does_span_coupling = (os%thermo_spans_coupling .and. (os%dt_therm > 1.5*dt_coupling))
593 
594  if (thermo_does_span_coupling) then
595  dt_therm = dt_coupling * floor(os%dt_therm / dt_coupling + 0.001)
596  nts = floor(dt_therm/dt_dyn + 0.001)
597  else
598  nts = max(1,min(n_max,floor(os%dt_therm/dt_dyn + 0.001)))
599  n_last_thermo = 0
600  endif
601 
602  time1 = time_seg_start ; t_elapsed_seg = 0.0
603  do n=1,n_max
604  if (os%diabatic_first) then
605  if (thermo_does_span_coupling) call mom_error(fatal, &
606  "MOM is not yet set up to have restarts that work with "//&
607  "THERMO_SPANS_COUPLING and DIABATIC_FIRST.")
608  if (modulo(n-1,nts)==0) then
609  dtdia = dt_dyn*min(nts,n_max-(n-1))
610  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dtdia, os%MOM_CSp, &
611  waves=os%Waves, do_dynamics=.false., do_thermodynamics=.true., &
612  start_cycle=(n==1), end_cycle=.false., cycle_length=dt_coupling)
613  endif
614 
615  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_dyn, os%MOM_CSp, &
616  waves=os%Waves, do_dynamics=.true., do_thermodynamics=.false., &
617  start_cycle=.false., end_cycle=(n==n_max), cycle_length=dt_coupling)
618  else
619  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_dyn, os%MOM_CSp, &
620  waves=os%Waves, do_dynamics=.true., do_thermodynamics=.false., &
621  start_cycle=(n==1), end_cycle=.false., cycle_length=dt_coupling)
622 
623  step_thermo = .false.
624  if (thermo_does_span_coupling) then
625  dtdia = dt_therm
626  step_thermo = mom_state_is_synchronized(os%MOM_CSp, adv_dyn=.true.)
627  elseif ((modulo(n,nts)==0) .or. (n==n_max)) then
628  dtdia = dt_dyn*(n - n_last_thermo)
629  n_last_thermo = n
630  step_thermo = .true.
631  endif
632 
633  if (step_thermo) then
634  ! Back up Time1 to the start of the thermodynamic segment.
635  time1 = time1 - real_to_time(dtdia - dt_dyn)
636  call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dtdia, os%MOM_CSp, &
637  waves=os%Waves, do_dynamics=.false., do_thermodynamics=.true., &
638  start_cycle=.false., end_cycle=(n==n_max), cycle_length=dt_coupling)
639  endif
640  endif
641 
642  t_elapsed_seg = t_elapsed_seg + dt_dyn
643  time1 = time_seg_start + real_to_time(t_elapsed_seg)
644  enddo
645  endif
646 
647  if (do_dyn) os%Time_dyn = time_seg_start + ocean_coupling_time_step
648  if (do_dyn) os%nstep = os%nstep + 1
649  os%Time = time_thermo_start ! Reset the clock to compensate for shared pointers.
650  if (do_thermo) os%Time = os%Time + ocean_coupling_time_step
651  if (do_thermo) os%nstep_thermo = os%nstep_thermo + 1
652 
653  if (do_dyn) then
654  call mech_forcing_diags(os%forces, dt_coupling, os%grid, os%Time_dyn, os%diag, os%forcing_CSp%handles)
655  endif
656 
657  if (os%fluxes%fluxes_used .and. do_thermo) then
658  call forcing_diagnostics(os%fluxes, os%sfc_state, os%grid, os%US, os%Time, os%diag, os%forcing_CSp%handles)
659  endif
660 
661 ! Translate state into Ocean.
662 ! call convert_state_to_ocean_type(OS%sfc_state, Ocean_sfc, OS%grid, OS%US, &
663 ! Ice_ocean_boundary%p, OS%press_to_z)
664  call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
665  time1 = os%Time ; if (do_dyn) time1 = os%Time_dyn
666  call coupler_type_send_data(ocean_sfc%fields, time1)
667 
668  call calltree_leave("update_ocean_model()")