MOM6
mom_surface_forcing Module Reference

Detailed Description

Functions that calculate the surface wind stresses and fluxes of buoyancy or temperature/salinity andfresh water, in ocean-only (solo) mode.

These functions are called every time step, even if the wind stresses or buoyancy fluxes are constant in time - in that case these routines return quickly without doing anything. In addition, any I/O of forcing fields is controlled by surface_forcing_init, located in this file.

Data Types

type  surface_forcing_cs
 Structure containing pointers to the forcing fields that may be used to drive MOM. All fluxes are positive into the ocean. More...
 
integer id_clock_forcing
 A CPU time clock. More...
 
subroutine, public set_forcing (sfc_state, forces, fluxes, day_start, day_interval, G, US, CS)
 Calls subroutines in this file to get surface forcing fields. More...
 
subroutine wind_forcing_const (sfc_state, forces, tau_x0, tau_y0, day, G, US, CS)
 Sets the surface wind stresses to constant values. More...
 
subroutine wind_forcing_2gyre (sfc_state, forces, day, G, US, CS)
 Sets the surface wind stresses to set up two idealized gyres. More...
 
subroutine wind_forcing_1gyre (sfc_state, forces, day, G, US, CS)
 Sets the surface wind stresses to set up a single idealized gyre. More...
 
subroutine wind_forcing_gyres (sfc_state, forces, day, G, US, CS)
 Sets the surface wind stresses to set up idealized gyres. More...
 
subroutine wind_forcing_from_file (sfc_state, forces, day, G, US, CS)
 
subroutine wind_forcing_by_data_override (sfc_state, forces, day, G, US, CS)
 
subroutine buoyancy_forcing_from_files (sfc_state, fluxes, day, dt, G, US, CS)
 Specifies zero surface bouyancy fluxes from input files. More...
 
subroutine buoyancy_forcing_from_data_override (sfc_state, fluxes, day, dt, G, US, CS)
 Specifies zero surface bouyancy fluxes from data over-ride. More...
 
subroutine buoyancy_forcing_zero (sfc_state, fluxes, day, dt, G, CS)
 This subroutine specifies zero surface bouyancy fluxes. More...
 
subroutine buoyancy_forcing_const (sfc_state, fluxes, day, dt, G, CS)
 Sets up spatially and temporally constant surface heat fluxes. More...
 
subroutine buoyancy_forcing_linear (sfc_state, fluxes, day, dt, G, US, CS)
 Sets surface fluxes of heat and salinity by restoring to temperature and salinity profiles that vary linearly with latitude. More...
 
subroutine, public forcing_save_restart (CS, G, Time, directory, time_stamped, filename_suffix)
 Save a restart file for the forcing fields. More...
 
subroutine, public surface_forcing_init (Time, G, US, param_file, diag, CS, tracer_flow_CSp)
 Initialize the surface forcing module. More...
 
subroutine surface_forcing_end (CS, fluxes)
 Deallocate memory associated with the surface forcing module. More...
 

Function/Subroutine Documentation

◆ buoyancy_forcing_const()

subroutine mom_surface_forcing::buoyancy_forcing_const ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(surface_forcing_cs), pointer  CS 
)
private

Sets up spatially and temporally constant surface heat fluxes.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]dayThe time of the fluxes
[in]dtThe amount of time over which the fluxes apply [s]
[in]gThe ocean's grid structure
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 1255 of file MOM_surface_forcing.F90.

1255  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1256  !! describe the surface state of the ocean.
1257  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1258  type(time_type), intent(in) :: day !< The time of the fluxes
1259  real, intent(in) :: dt !< The amount of time over which
1260  !! the fluxes apply [s]
1261  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1262  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1263  !! a previous surface_forcing_init call
1264  ! Local variables
1265  integer :: i, j, is, ie, js, je
1266  call calltree_enter("buoyancy_forcing_const, MOM_surface_forcing.F90")
1267  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1268 
1269  if (cs%use_temperature) then
1270  do j=js,je ; do i=is,ie
1271  fluxes%evap(i,j) = 0.0
1272  fluxes%lprec(i,j) = 0.0
1273  fluxes%fprec(i,j) = 0.0
1274  fluxes%vprec(i,j) = 0.0
1275  fluxes%lrunoff(i,j) = 0.0
1276  fluxes%frunoff(i,j) = 0.0
1277  fluxes%lw(i,j) = 0.0
1278  fluxes%latent(i,j) = 0.0
1279  fluxes%sens(i,j) = cs%constantHeatForcing * g%mask2dT(i,j)
1280  fluxes%sw(i,j) = 0.0
1281  fluxes%latent_evap_diag(i,j) = 0.0
1282  fluxes%latent_fprec_diag(i,j) = 0.0
1283  fluxes%latent_frunoff_diag(i,j) = 0.0
1284  enddo ; enddo
1285  else
1286  do j=js,je ; do i=is,ie
1287  fluxes%buoy(i,j) = 0.0
1288  enddo ; enddo
1289  endif
1290 
1291  call calltree_leave("buoyancy_forcing_const")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ buoyancy_forcing_from_data_override()

subroutine mom_surface_forcing::buoyancy_forcing_from_data_override ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Specifies zero surface bouyancy fluxes from data over-ride.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]dayThe time of the fluxes
[in]dtThe amount of time over which the fluxes apply [s]
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 1035 of file MOM_surface_forcing.F90.

1035  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1036  !! describe the surface state of the ocean.
1037  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1038  type(time_type), intent(in) :: day !< The time of the fluxes
1039  real, intent(in) :: dt !< The amount of time over which
1040  !! the fluxes apply [s]
1041  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1042  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1043  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1044  !! a previous surface_forcing_init call
1045  ! Local variables
1046  real, dimension(SZI_(G),SZJ_(G)) :: &
1047  temp, & ! A 2-d temporary work array with various units.
1048  SST_anom, & ! Instantaneous sea surface temperature anomalies from a
1049  ! target (observed) value [degC].
1050  sss_anom, & ! Instantaneous sea surface salinity anomalies from a target
1051  ! (observed) value [ppt].
1052  sss_mean ! A (mean?) salinity about which to normalize local salinity
1053  ! anomalies when calculating restorative precipitation
1054  ! anomalies [ppt].
1055  real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
1056  ! mass fluxes [R Z s m2 kg-1 T-1 ~> 1].
1057  real :: rhoXcp ! The mean density times the heat capacity [J m-3 degC-1].
1058 
1059  integer :: time_lev_daily ! The time levels to read for fields with
1060  integer :: time_lev_monthly ! daily and montly cycles.
1061  integer :: itime_lev ! The time level that is used for a field.
1062 
1063  integer :: days, seconds
1064  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
1065  integer :: is_in, ie_in, js_in, je_in
1066 
1067  call calltree_enter("buoyancy_forcing_from_data_override, MOM_surface_forcing.F90")
1068 
1069  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1070  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1071  kg_m2_s_conversion = us%kg_m3_to_R*us%m_to_Z*us%T_to_s
1072 
1073  if (cs%use_temperature) rhoxcp = cs%Rho0 * fluxes%C_p
1074 
1075  if (.not.cs%dataOverrideIsInitialized) then
1076  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
1077  cs%dataOverrideIsInitialized = .true.
1078  endif
1079 
1080  is_in = g%isc - g%isd + 1
1081  ie_in = g%iec - g%isd + 1
1082  js_in = g%jsc - g%jsd + 1
1083  je_in = g%jec - g%jsd + 1
1084 
1085  call data_override('OCN', 'lw', fluxes%LW(:,:), day, &
1086  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1087  call data_override('OCN', 'evap', fluxes%evap(:,:), day, &
1088  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1089 
1090  ! note the sign convention
1091  do j=js,je ; do i=is,ie
1092  ! This is dangerous because it is not clear whether the data files have been read!
1093  fluxes%evap(i,j) = -fluxes%evap(i,j) ! Normal convention is positive into the ocean
1094  ! but evap is normally a positive quantity in the files
1095  fluxes%latent(i,j) = cs%latent_heat_vapor*fluxes%evap(i,j)
1096  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
1097  fluxes%evap(i,j) = kg_m2_s_conversion*fluxes%evap(i,j)
1098  enddo ; enddo
1099 
1100  call data_override('OCN', 'sens', fluxes%sens(:,:), day, &
1101  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1102 
1103  ! note the sign convention
1104  do j=js,je ; do i=is,ie
1105  fluxes%sens(i,j) = -fluxes%sens(i,j) ! Normal convention is positive into the ocean
1106  ! but sensible is normally a positive quantity in the files
1107  enddo ; enddo
1108 
1109  call data_override('OCN', 'sw', fluxes%sw(:,:), day, &
1110  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1111 
1112  call data_override('OCN', 'snow', fluxes%fprec(:,:), day, &
1113  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1114 
1115  call data_override('OCN', 'rain', fluxes%lprec(:,:), day, &
1116  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1117 
1118  call data_override('OCN', 'runoff', fluxes%lrunoff(:,:), day, &
1119  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1120 
1121  call data_override('OCN', 'calving', fluxes%frunoff(:,:), day, &
1122  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in) ! scale=kg_m2_s_conversion
1123 
1124  if (kg_m2_s_conversion /= 1.0) then ; do j=js,je ; do i=is,ie
1125  fluxes%lprec(i,j) = fluxes%lprec(i,j) * kg_m2_s_conversion
1126  fluxes%fprec(i,j) = fluxes%fprec(i,j) * kg_m2_s_conversion
1127  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * kg_m2_s_conversion
1128  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * kg_m2_s_conversion
1129  enddo ; enddo ; endif
1130 
1131 ! Read the SST and SSS fields for damping.
1132  if (cs%restorebuoy) then !#CTRL# .or. associated(CS%ctrl_forcing_CSp)) then
1133  call data_override('OCN', 'SST_restore', cs%T_restore(:,:), day, &
1134  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1135 
1136  call data_override('OCN', 'SSS_restore', cs%S_restore(:,:), day, &
1137  is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
1138 
1139  endif
1140 
1141  ! restoring boundary fluxes
1142  if (cs%restorebuoy) then
1143  if (cs%use_temperature) then
1144  do j=js,je ; do i=is,ie
1145  if (g%mask2dT(i,j) > 0) then
1146  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
1147  ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
1148  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
1149  (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
1150  (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
1151  else
1152  fluxes%heat_added(i,j) = 0.0
1153  fluxes%vprec(i,j) = 0.0
1154  endif
1155  enddo ; enddo
1156  else
1157  do j=js,je ; do i=is,ie
1158  if (g%mask2dT(i,j) > 0) then
1159  fluxes%buoy(i,j) = us%kg_m3_to_R * (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1160  (cs%G_Earth * cs%Flux_const / cs%Rho0)
1161  else
1162  fluxes%buoy(i,j) = 0.0
1163  endif
1164  enddo ; enddo
1165  endif
1166  else ! not RESTOREBUOY
1167  if (.not.cs%use_temperature) then
1168  call mom_error(fatal, "buoyancy_forcing in MOM_surface_forcing: "// &
1169  "The fluxes need to be defined without RESTOREBUOY.")
1170  endif
1171  endif ! end RESTOREBUOY
1172 
1173 
1174  ! mask out land points and compute heat content of water fluxes
1175  ! assume liquid precip enters ocean at SST
1176  ! assume frozen precip enters ocean at 0degC
1177  ! assume liquid runoff enters ocean at SST
1178  ! assume solid runoff (calving) enters ocean at 0degC
1179  ! mass leaving ocean has heat_content determined in MOM_diabatic_driver.F90
1180  do j=js,je ; do i=is,ie
1181  fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
1182  fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
1183  fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
1184  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
1185  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
1186  fluxes%LW(i,j) = fluxes%LW(i,j) * g%mask2dT(i,j)
1187  fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
1188  fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
1189  fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
1190 
1191  fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
1192  fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
1193  fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
1194  enddo ; enddo
1195 
1196 
1197 !#CTRL# if (associated(CS%ctrl_forcing_CSp)) then
1198 !#CTRL# do j=js,je ; do i=is,ie
1199 !#CTRL# SST_anom(i,j) = sfc_state%SST(i,j) - CS%T_Restore(i,j)
1200 !#CTRL# SSS_anom(i,j) = sfc_state%SSS(i,j) - CS%S_Restore(i,j)
1201 !#CTRL# SSS_mean(i,j) = 0.5*(sfc_state%SSS(i,j) + CS%S_Restore(i,j))
1202 !#CTRL# enddo ; enddo
1203 !#CTRL# call apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, fluxes%heat_added, &
1204 !#CTRL# fluxes%vprec, day, dt, G, CS%ctrl_forcing_CSp)
1205 !#CTRL# endif
1206 
1207  call calltree_leave("buoyancy_forcing_from_data_override")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ buoyancy_forcing_from_files()

subroutine mom_surface_forcing::buoyancy_forcing_from_files ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Specifies zero surface bouyancy fluxes from input files.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]dayThe time of the fluxes
[in]dtThe amount of time over which the fluxes apply [s]
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 755 of file MOM_surface_forcing.F90.

755  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
756  !! describe the surface state of the ocean.
757  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
758  type(time_type), intent(in) :: day !< The time of the fluxes
759  real, intent(in) :: dt !< The amount of time over which
760  !! the fluxes apply [s]
761  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
762  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
763  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
764  !! a previous surface_forcing_init call
765  ! Local variables
766  real, dimension(SZI_(G),SZJ_(G)) :: &
767  temp, & ! A 2-d temporary work array with various units.
768  SST_anom, & ! Instantaneous sea surface temperature anomalies from a
769  ! target (observed) value [degC].
770  sss_anom, & ! Instantaneous sea surface salinity anomalies from a target
771  ! (observed) value [ppt].
772  sss_mean ! A (mean?) salinity about which to normalize local salinity
773  ! anomalies when calculating restorative precipitation
774  ! anomalies [ppt].
775 
776  real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
777  ! mass fluxes [R Z s m2 kg-1 T-1 ~> 1].
778  real :: rhoXcp ! reference density times heat capacity [J m-3 degC-1]
779 
780  integer :: time_lev_daily ! time levels to read for fields with daily cycle
781  integer :: time_lev_monthly ! time levels to read for fields with monthly cycle
782  integer :: time_lev ! time level that for a field
783 
784  integer :: days, seconds
785  integer :: i, j, is, ie, js, je
786 
787  call calltree_enter("buoyancy_forcing_from_files, MOM_surface_forcing.F90")
788 
789  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
790  kg_m2_s_conversion = us%kg_m3_to_R*us%m_to_Z*us%T_to_s
791 
792  if (cs%use_temperature) rhoxcp = us%R_to_kg_m3*cs%Rho0 * fluxes%C_p
793 
794  ! Read the buoyancy forcing file
795  call get_time(day, seconds, days)
796 
797  time_lev_daily = days - 365*floor(real(days) / 365.0)
798 
799  if (time_lev_daily < 31) then ; time_lev_monthly = 0
800  elseif (time_lev_daily < 59) then ; time_lev_monthly = 1
801  elseif (time_lev_daily < 90) then ; time_lev_monthly = 2
802  elseif (time_lev_daily < 120) then ; time_lev_monthly = 3
803  elseif (time_lev_daily < 151) then ; time_lev_monthly = 4
804  elseif (time_lev_daily < 181) then ; time_lev_monthly = 5
805  elseif (time_lev_daily < 212) then ; time_lev_monthly = 6
806  elseif (time_lev_daily < 243) then ; time_lev_monthly = 7
807  elseif (time_lev_daily < 273) then ; time_lev_monthly = 8
808  elseif (time_lev_daily < 304) then ; time_lev_monthly = 9
809  elseif (time_lev_daily < 334) then ; time_lev_monthly = 10
810  else ; time_lev_monthly = 11
811  endif
812 
813  time_lev_daily = time_lev_daily +1
814  time_lev_monthly = time_lev_monthly+1
815 
816  if (time_lev_daily /= cs%buoy_last_lev_read) then
817 
818  ! longwave
819  select case (cs%LW_nlev)
820  case (12) ; time_lev = time_lev_monthly
821  case (365) ; time_lev = time_lev_daily
822  case default ; time_lev = 1
823  end select
824  call mom_read_data(cs%longwave_file, cs%LW_var, fluxes%LW(:,:), &
825  g%Domain, timelevel=time_lev)
826  if (cs%archaic_OMIP_file) then
827  call mom_read_data(cs%longwaveup_file, "lwup_sfc", temp(:,:), g%Domain, &
828  timelevel=time_lev)
829  do j=js,je ; do i=is,ie ; fluxes%LW(i,j) = fluxes%LW(i,j) - temp(i,j) ; enddo ; enddo
830  endif
831  cs%LW_last_lev = time_lev
832 
833  ! evaporation
834  select case (cs%evap_nlev)
835  case (12) ; time_lev = time_lev_monthly
836  case (365) ; time_lev = time_lev_daily
837  case default ; time_lev = 1
838  end select
839  if (cs%archaic_OMIP_file) then
840  call mom_read_data(cs%evaporation_file, cs%evap_var, temp(:,:), &
841  g%Domain, timelevel=time_lev)
842  do j=js,je ; do i=is,ie
843  fluxes%latent(i,j) = -cs%latent_heat_vapor*temp(i,j)
844  fluxes%evap(i,j) = -kg_m2_s_conversion*temp(i,j)
845  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
846  enddo ; enddo
847  else
848  call mom_read_data(cs%evaporation_file, cs%evap_var, fluxes%evap(:,:), &
849  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
850  endif
851  cs%evap_last_lev = time_lev
852 
853  select case (cs%latent_nlev)
854  case (12) ; time_lev = time_lev_monthly
855  case (365) ; time_lev = time_lev_daily
856  case default ; time_lev = 1
857  end select
858  if (.not.cs%archaic_OMIP_file) then
859  call mom_read_data(cs%latentheat_file, cs%latent_var, fluxes%latent(:,:), &
860  g%Domain, timelevel=time_lev)
861  do j=js,je ; do i=is,ie
862  fluxes%latent_evap_diag(i,j) = fluxes%latent(i,j)
863  enddo ; enddo
864  endif
865  cs%latent_last_lev = time_lev
866 
867  select case (cs%sens_nlev)
868  case (12) ; time_lev = time_lev_monthly
869  case (365) ; time_lev = time_lev_daily
870  case default ; time_lev = 1
871  end select
872  if (cs%archaic_OMIP_file) then
873  call mom_read_data(cs%sensibleheat_file, cs%sens_var, temp(:,:), &
874  g%Domain, timelevel=time_lev)
875  do j=js,je ; do i=is,ie ; fluxes%sens(i,j) = -temp(i,j) ; enddo ; enddo
876  else
877  call mom_read_data(cs%sensibleheat_file, cs%sens_var, fluxes%sens(:,:), &
878  g%Domain, timelevel=time_lev)
879  endif
880  cs%sens_last_lev = time_lev
881 
882  select case (cs%SW_nlev)
883  case (12) ; time_lev = time_lev_monthly
884  case (365) ; time_lev = time_lev_daily
885  case default ; time_lev = 1
886  end select
887  call mom_read_data(cs%shortwave_file, cs%SW_var, fluxes%sw(:,:), &
888  g%Domain, timelevel=time_lev)
889  if (cs%archaic_OMIP_file) then
890  call mom_read_data(cs%shortwaveup_file, "swup_sfc", temp(:,:), &
891  g%Domain, timelevel=time_lev)
892  do j=js,je ; do i=is,ie
893  fluxes%sw(i,j) = fluxes%sw(i,j) - temp(i,j)
894  enddo ; enddo
895  endif
896  cs%SW_last_lev = time_lev
897 
898  select case (cs%precip_nlev)
899  case (12) ; time_lev = time_lev_monthly
900  case (365) ; time_lev = time_lev_daily
901  case default ; time_lev = 1
902  end select
903  call mom_read_data(cs%snow_file, cs%snow_var, &
904  fluxes%fprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
905  call mom_read_data(cs%rain_file, cs%rain_var, &
906  fluxes%lprec(:,:), g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
907  if (cs%archaic_OMIP_file) then
908  do j=js,je ; do i=is,ie
909  fluxes%lprec(i,j) = fluxes%lprec(i,j) - fluxes%fprec(i,j)
910  enddo ; enddo
911  endif
912  cs%precip_last_lev = time_lev
913 
914  select case (cs%runoff_nlev)
915  case (12) ; time_lev = time_lev_monthly
916  case (365) ; time_lev = time_lev_daily
917  case default ; time_lev = 1
918  end select
919  if (cs%archaic_OMIP_file) then
920  call mom_read_data(cs%runoff_file, cs%lrunoff_var, temp(:,:), &
921  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
922  do j=js,je ; do i=is,ie
923  fluxes%lrunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
924  enddo ; enddo
925  call mom_read_data(cs%runoff_file, cs%frunoff_var, temp(:,:), &
926  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
927  do j=js,je ; do i=is,ie
928  fluxes%frunoff(i,j) = temp(i,j)*us%m_to_L**2*g%IareaT(i,j)
929  enddo ; enddo
930  else
931  call mom_read_data(cs%runoff_file, cs%lrunoff_var, fluxes%lrunoff(:,:), &
932  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
933  call mom_read_data(cs%runoff_file, cs%frunoff_var, fluxes%frunoff(:,:), &
934  g%Domain, timelevel=time_lev, scale=kg_m2_s_conversion)
935  endif
936  cs%runoff_last_lev = time_lev
937 
938 ! Read the SST and SSS fields for damping.
939  if (cs%restorebuoy) then !#CTRL# .or. associated(CS%ctrl_forcing_CSp)) then
940  select case (cs%SST_nlev)
941  case (12) ; time_lev = time_lev_monthly
942  case (365) ; time_lev = time_lev_daily
943  case default ; time_lev = 1
944  end select
945  call mom_read_data(cs%SSTrestore_file, cs%SST_restore_var, &
946  cs%T_Restore(:,:), g%Domain, timelevel=time_lev)
947  cs%SST_last_lev = time_lev
948 
949  select case (cs%SSS_nlev)
950  case (12) ; time_lev = time_lev_monthly
951  case (365) ; time_lev = time_lev_daily
952  case default ; time_lev = 1
953  end select
954  call mom_read_data(cs%salinityrestore_file, cs%SSS_restore_var, &
955  cs%S_Restore(:,:), g%Domain, timelevel=time_lev)
956  cs%SSS_last_lev = time_lev
957  endif
958  cs%buoy_last_lev_read = time_lev_daily
959 
960  ! mask out land points and compute heat content of water fluxes
961  ! assume liquid precip enters ocean at SST
962  ! assume frozen precip enters ocean at 0degC
963  ! assume liquid runoff enters ocean at SST
964  ! assume solid runoff (calving) enters ocean at 0degC
965  ! mass leaving the ocean has heat_content determined in MOM_diabatic_driver.F90
966  do j=js,je ; do i=is,ie
967  fluxes%evap(i,j) = fluxes%evap(i,j) * g%mask2dT(i,j)
968  fluxes%lprec(i,j) = fluxes%lprec(i,j) * g%mask2dT(i,j)
969  fluxes%fprec(i,j) = fluxes%fprec(i,j) * g%mask2dT(i,j)
970  fluxes%lrunoff(i,j) = fluxes%lrunoff(i,j) * g%mask2dT(i,j)
971  fluxes%frunoff(i,j) = fluxes%frunoff(i,j) * g%mask2dT(i,j)
972  fluxes%LW(i,j) = fluxes%LW(i,j) * g%mask2dT(i,j)
973  fluxes%sens(i,j) = fluxes%sens(i,j) * g%mask2dT(i,j)
974  fluxes%sw(i,j) = fluxes%sw(i,j) * g%mask2dT(i,j)
975  fluxes%latent(i,j) = fluxes%latent(i,j) * g%mask2dT(i,j)
976 
977  fluxes%latent_evap_diag(i,j) = fluxes%latent_evap_diag(i,j) * g%mask2dT(i,j)
978  fluxes%latent_fprec_diag(i,j) = -fluxes%fprec(i,j)*cs%latent_heat_fusion
979  fluxes%latent_frunoff_diag(i,j) = -fluxes%frunoff(i,j)*cs%latent_heat_fusion
980  enddo ; enddo
981 
982  endif ! time_lev /= CS%buoy_last_lev_read
983 
984 
985  ! restoring surface boundary fluxes
986  if (cs%restorebuoy) then
987 
988  if (cs%use_temperature) then
989  do j=js,je ; do i=is,ie
990  if (g%mask2dT(i,j) > 0) then
991  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
992  ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const_T)
993  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const_S) * &
994  (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
995  (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
996  else
997  fluxes%heat_added(i,j) = 0.0
998  fluxes%vprec(i,j) = 0.0
999  endif
1000  enddo ; enddo
1001  else
1002  do j=js,je ; do i=is,ie
1003  if (g%mask2dT(i,j) > 0) then
1004  fluxes%buoy(i,j) = us%kg_m3_to_R * (cs%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1005  (cs%G_Earth * cs%Flux_const / cs%Rho0)
1006  else
1007  fluxes%buoy(i,j) = 0.0
1008  endif
1009  enddo ; enddo
1010  endif
1011 
1012  else ! not RESTOREBUOY
1013  if (.not.cs%use_temperature) then
1014  call mom_error(fatal, "buoyancy_forcing in MOM_surface_forcing: "// &
1015  "The fluxes need to be defined without RESTOREBUOY.")
1016  endif
1017 
1018  endif ! end RESTOREBUOY
1019 
1020 !#CTRL# if (associated(CS%ctrl_forcing_CSp)) then
1021 !#CTRL# do j=js,je ; do i=is,ie
1022 !#CTRL# SST_anom(i,j) = sfc_state%SST(i,j) - CS%T_Restore(i,j)
1023 !#CTRL# SSS_anom(i,j) = sfc_state%SSS(i,j) - CS%S_Restore(i,j)
1024 !#CTRL# SSS_mean(i,j) = 0.5*(sfc_state%SSS(i,j) + CS%S_Restore(i,j))
1025 !#CTRL# enddo ; enddo
1026 !#CTRL# call apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, fluxes%heat_added, &
1027 !#CTRL# fluxes%vprec, day, dt, G, CS%ctrl_forcing_CSp)
1028 !#CTRL# endif
1029 
1030  call calltree_leave("buoyancy_forcing_from_files")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ buoyancy_forcing_linear()

subroutine mom_surface_forcing::buoyancy_forcing_linear ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Sets surface fluxes of heat and salinity by restoring to temperature and salinity profiles that vary linearly with latitude.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]dayThe time of the fluxes
[in]dtThe amount of time over which the fluxes apply [s]
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 1297 of file MOM_surface_forcing.F90.

1297  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1298  !! describe the surface state of the ocean.
1299  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1300  type(time_type), intent(in) :: day !< The time of the fluxes
1301  real, intent(in) :: dt !< The amount of time over which
1302  !! the fluxes apply [s]
1303  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1304  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1305  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1306  !! a previous surface_forcing_init call
1307  ! Local variables
1308  real :: y, T_restore, S_restore
1309  integer :: i, j, is, ie, js, je
1310 
1311  call calltree_enter("buoyancy_forcing_linear, MOM_surface_forcing.F90")
1312  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1313 
1314  ! This case has no surface buoyancy forcing.
1315  if (cs%use_temperature) then
1316  do j=js,je ; do i=is,ie
1317  fluxes%evap(i,j) = 0.0
1318  fluxes%lprec(i,j) = 0.0
1319  fluxes%fprec(i,j) = 0.0
1320  fluxes%vprec(i,j) = 0.0
1321  fluxes%lrunoff(i,j) = 0.0
1322  fluxes%frunoff(i,j) = 0.0
1323  fluxes%lw(i,j) = 0.0
1324  fluxes%latent(i,j) = 0.0
1325  fluxes%sens(i,j) = 0.0
1326  fluxes%sw(i,j) = 0.0
1327  fluxes%latent_evap_diag(i,j) = 0.0
1328  fluxes%latent_fprec_diag(i,j) = 0.0
1329  fluxes%latent_frunoff_diag(i,j) = 0.0
1330  enddo ; enddo
1331  else
1332  do j=js,je ; do i=is,ie
1333  fluxes%buoy(i,j) = 0.0
1334  enddo ; enddo
1335  endif
1336 
1337  if (cs%restorebuoy) then
1338  if (cs%use_temperature) then
1339  do j=js,je ; do i=is,ie
1340  y = (g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat
1341  t_restore = cs%T_south + (cs%T_north-cs%T_south)*y
1342  s_restore = cs%S_south + (cs%S_north-cs%S_south)*y
1343  if (g%mask2dT(i,j) > 0) then
1344  fluxes%heat_added(i,j) = g%mask2dT(i,j) * (us%R_to_kg_m3*us%Z_to_m*us%s_to_T) * &
1345  ((t_restore - sfc_state%SST(i,j)) * ((cs%Rho0 * fluxes%C_p) * cs%Flux_const))
1346  fluxes%vprec(i,j) = - (cs%Rho0*cs%Flux_const) * &
1347  (s_restore - sfc_state%SSS(i,j)) / &
1348  (0.5*(sfc_state%SSS(i,j) + s_restore))
1349  else
1350  fluxes%heat_added(i,j) = 0.0
1351  fluxes%vprec(i,j) = 0.0
1352  endif
1353  enddo ; enddo
1354  else
1355  call mom_error(fatal, "buoyancy_forcing_linear in MOM_surface_forcing: "// &
1356  "RESTOREBUOY to linear not written yet.")
1357  !do j=js,je ; do i=is,ie
1358  ! if (G%mask2dT(i,j) > 0) then
1359  ! fluxes%buoy(i,j) = US%kg_m3_to_R * (CS%Dens_Restore(i,j) - sfc_state%sfc_density(i,j)) * &
1360  ! (CS%G_Earth * CS%Flux_const / CS%Rho0)
1361  ! else
1362  ! fluxes%buoy(i,j) = 0.0
1363  ! endif
1364  !enddo ; enddo
1365  endif
1366  else ! not RESTOREBUOY
1367  if (.not.cs%use_temperature) then
1368  call mom_error(fatal, "buoyancy_forcing_linear in MOM_surface_forcing: "// &
1369  "The fluxes need to be defined without RESTOREBUOY.")
1370  endif
1371  endif ! end RESTOREBUOY
1372 
1373  call calltree_leave("buoyancy_forcing_linear")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ buoyancy_forcing_zero()

subroutine mom_surface_forcing::buoyancy_forcing_zero ( type(surface), intent(inout)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(surface_forcing_cs), pointer  CS 
)
private

This subroutine specifies zero surface bouyancy fluxes.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]dayThe time of the fluxes
[in]dtThe amount of time over which the fluxes apply [s]
[in]gThe ocean's grid structure
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 1212 of file MOM_surface_forcing.F90.

1212  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
1213  !! describe the surface state of the ocean.
1214  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1215  type(time_type), intent(in) :: day !< The time of the fluxes
1216  real, intent(in) :: dt !< The amount of time over which
1217  !! the fluxes apply [s]
1218  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1219  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1220  !! a previous surface_forcing_init call
1221  ! Local variables
1222  integer :: i, j, is, ie, js, je
1223 
1224  call calltree_enter("buoyancy_forcing_zero, MOM_surface_forcing.F90")
1225  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1226 
1227  if (cs%use_temperature) then
1228  do j=js,je ; do i=is,ie
1229  fluxes%evap(i,j) = 0.0
1230  fluxes%lprec(i,j) = 0.0
1231  fluxes%fprec(i,j) = 0.0
1232  fluxes%vprec(i,j) = 0.0
1233  fluxes%lrunoff(i,j) = 0.0
1234  fluxes%frunoff(i,j) = 0.0
1235  fluxes%lw(i,j) = 0.0
1236  fluxes%latent(i,j) = 0.0
1237  fluxes%sens(i,j) = 0.0
1238  fluxes%sw(i,j) = 0.0
1239  fluxes%latent_evap_diag(i,j) = 0.0
1240  fluxes%latent_fprec_diag(i,j) = 0.0
1241  fluxes%latent_frunoff_diag(i,j) = 0.0
1242  enddo ; enddo
1243  else
1244  do j=js,je ; do i=is,ie
1245  fluxes%buoy(i,j) = 0.0
1246  enddo ; enddo
1247  endif
1248 
1249  call calltree_leave("buoyancy_forcing_zero")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ forcing_save_restart()

subroutine, public mom_surface_forcing::forcing_save_restart ( type(surface_forcing_cs), pointer  CS,
type(ocean_grid_type), intent(inout)  G,
type(time_type), intent(in)  Time,
character(len=*), intent(in)  directory,
logical, intent(in), optional  time_stamped,
character(len=*), intent(in), optional  filename_suffix 
)

Save a restart file for the forcing fields.

Parameters
cspointer to control struct returned by a previous surface_forcing_init call
[in,out]gThe ocean's grid structure
[in]timemodel time at this call; needed for mpp_write calls
[in]directorydirectory into which to write these restart files
[in]time_stampedIf true, the restart file names include a unique time stamp; the default is false.
[in]filename_suffixoptional suffix (e.g., a time-stamp) to append to the restart fname

Definition at line 1379 of file MOM_surface_forcing.F90.

1379  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1380  !! a previous surface_forcing_init call
1381  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1382  type(time_type), intent(in) :: Time !< model time at this call; needed for mpp_write calls
1383  character(len=*), intent(in) :: directory !< directory into which to write these restart files
1384  logical, optional, intent(in) :: time_stamped !< If true, the restart file names
1385  !! include a unique time stamp; the default is false.
1386  character(len=*), optional, intent(in) :: filename_suffix !< optional suffix (e.g., a time-stamp)
1387  !! to append to the restart fname
1388 
1389  if (.not.associated(cs)) return
1390  if (.not.associated(cs%restart_CSp)) return
1391 
1392  call save_restart(directory, time, g, cs%restart_CSp, time_stamped)
1393 

References mom_restart::save_restart().

Referenced by mom_main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_forcing()

subroutine, public mom_surface_forcing::set_forcing ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day_start,
type(time_type), intent(in)  day_interval,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)

Calls subroutines in this file to get surface forcing fields.

It also allocates and initializes the fields in the forcing and mech_forcing types the first time it is called.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in,out]fluxesA structure containing thermodynamic forcing fields
[in]day_startThe start time of the fluxes
[in]day_intervalLength of time over which these fluxes applied
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 221 of file MOM_surface_forcing.F90.

221  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
222  !! describe the surface state of the ocean.
223  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
224  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
225  type(time_type), intent(in) :: day_start !< The start time of the fluxes
226  type(time_type), intent(in) :: day_interval !< Length of time over which these fluxes applied
227  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
228  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
229  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
230  !! a previous surface_forcing_init call
231  ! Local variables
232  real :: dt ! length of time over which fluxes applied [s]
233  type(time_type) :: day_center ! central time of the fluxes.
234  integer :: isd, ied, jsd, jed
235  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
236 
237  call cpu_clock_begin(id_clock_forcing)
238  call calltree_enter("set_forcing, MOM_surface_forcing.F90")
239 
240  day_center = day_start + day_interval/2
241  dt = time_type_to_real(day_interval)
242 
243  if (cs%first_call_set_forcing) then
244  ! Allocate memory for the mechanical and thermodyanmic forcing fields.
245  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true., press=.true.)
246 
247  call allocate_forcing_type(g, fluxes, ustar=.true.)
248  if (trim(cs%buoy_config) /= "NONE") then
249  if ( cs%use_temperature ) then
250  call allocate_forcing_type(g, fluxes, water=.true., heat=.true., press=.true.)
251  if (cs%restorebuoy) then
252  call safe_alloc_ptr(cs%T_Restore,isd, ied, jsd, jed)
253  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
254  call safe_alloc_ptr(cs%S_Restore, isd, ied, jsd, jed)
255  endif
256  else ! CS%use_temperature false.
257  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
258 
259  if (cs%restorebuoy) call safe_alloc_ptr(cs%Dens_Restore, isd, ied, jsd, jed)
260  endif ! endif for CS%use_temperature
261  endif
262  endif
263 
264  ! calls to various wind options
265  if (cs%variable_winds .or. cs%first_call_set_forcing) then
266 
267  if (trim(cs%wind_config) == "file") then
268  call wind_forcing_from_file(sfc_state, forces, day_center, g, us, cs)
269  elseif (trim(cs%wind_config) == "data_override") then
270  call wind_forcing_by_data_override(sfc_state, forces, day_center, g, us, cs)
271  elseif (trim(cs%wind_config) == "2gyre") then
272  call wind_forcing_2gyre(sfc_state, forces, day_center, g, us, cs)
273  elseif (trim(cs%wind_config) == "1gyre") then
274  call wind_forcing_1gyre(sfc_state, forces, day_center, g, us, cs)
275  elseif (trim(cs%wind_config) == "gyres") then
276  call wind_forcing_gyres(sfc_state, forces, day_center, g, us, cs)
277  elseif (trim(cs%wind_config) == "zero") then
278  call wind_forcing_const(sfc_state, forces, 0., 0., day_center, g, us, cs)
279  elseif (trim(cs%wind_config) == "const") then
280  call wind_forcing_const(sfc_state, forces, cs%tau_x0, cs%tau_y0, day_center, g, us, cs)
281  elseif (trim(cs%wind_config) == "Neverland") then
282  call neverland_wind_forcing(sfc_state, forces, day_center, g, us, cs%Neverland_forcing_CSp)
283  elseif (trim(cs%wind_config) == "ideal_hurr") then
284  call idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
285  elseif (trim(cs%wind_config) == "SCM_ideal_hurr") then
286  call scm_idealized_hurricane_wind_forcing(sfc_state, forces, day_center, g, us, cs%idealized_hurricane_CSp)
287  elseif (trim(cs%wind_config) == "SCM_CVmix_tests") then
288  call scm_cvmix_tests_wind_forcing(sfc_state, forces, day_center, g, us, cs%SCM_CVmix_tests_CSp)
289  elseif (trim(cs%wind_config) == "USER") then
290  call user_wind_forcing(sfc_state, forces, day_center, g, us, cs%user_forcing_CSp)
291  elseif (cs%variable_winds .and. .not.cs%first_call_set_forcing) then
292  call mom_error(fatal, &
293  "MOM_surface_forcing: Variable winds defined with no wind config")
294  else
295  call mom_error(fatal, &
296  "MOM_surface_forcing:Unrecognized wind config "//trim(cs%wind_config))
297  endif
298  endif
299 
300  ! calls to various buoyancy forcing options
301  if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
302  (.not.cs%adiabatic)) then
303  if (trim(cs%buoy_config) == "file") then
304  call buoyancy_forcing_from_files(sfc_state, fluxes, day_center, dt, g, us, cs)
305  elseif (trim(cs%buoy_config) == "data_override") then
306  call buoyancy_forcing_from_data_override(sfc_state, fluxes, day_center, dt, g, us, cs)
307  elseif (trim(cs%buoy_config) == "zero") then
308  call buoyancy_forcing_zero(sfc_state, fluxes, day_center, dt, g, cs)
309  elseif (trim(cs%buoy_config) == "const") then
310  call buoyancy_forcing_const(sfc_state, fluxes, day_center, dt, g, cs)
311  elseif (trim(cs%buoy_config) == "linear") then
312  call buoyancy_forcing_linear(sfc_state, fluxes, day_center, dt, g, us, cs)
313  elseif (trim(cs%buoy_config) == "MESO") then
314  call meso_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%MESO_forcing_CSp)
315  elseif (trim(cs%buoy_config) == "Neverland") then
316  call neverland_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%Neverland_forcing_CSp)
317  elseif (trim(cs%buoy_config) == "SCM_CVmix_tests") then
318  call scm_cvmix_tests_buoyancy_forcing(sfc_state, fluxes, day_center, g, us, cs%SCM_CVmix_tests_CSp)
319  elseif (trim(cs%buoy_config) == "USER") then
320  call user_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%user_forcing_CSp)
321  elseif (trim(cs%buoy_config) == "BFB") then
322  call bfb_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%BFB_forcing_CSp)
323  elseif (trim(cs%buoy_config) == "dumbbell") then
324  call dumbbell_buoyancy_forcing(sfc_state, fluxes, day_center, dt, g, us, cs%dumbbell_forcing_CSp)
325  elseif (trim(cs%buoy_config) == "NONE") then
326  call mom_mesg("MOM_surface_forcing: buoyancy forcing has been set to omitted.")
327  elseif (cs%variable_buoyforce .and. .not.cs%first_call_set_forcing) then
328  call mom_error(fatal, &
329  "MOM_surface_forcing: Variable buoy defined with no buoy config.")
330  else
331  call mom_error(fatal, &
332  "MOM_surface_forcing: Unrecognized buoy config "//trim(cs%buoy_config))
333  endif
334  endif
335 
336  if (associated(cs%tracer_flow_CSp)) then
337  call call_tracer_set_forcing(sfc_state, fluxes, day_start, day_interval, g, cs%tracer_flow_CSp)
338  endif
339 
340  ! Allow for user-written code to alter the fluxes after all the above
341  call user_alter_forcing(sfc_state, fluxes, day_center, g, cs%urf_CS)
342 
343  ! Fields that exist in both the forcing and mech_forcing types must be copied.
344  if (cs%variable_winds .or. cs%first_call_set_forcing) then
345  call copy_common_forcing_fields(forces, fluxes, g)
346  call set_derived_forcing_fields(forces, fluxes, g, us, cs%Rho0)
347  endif
348 
349  if ((cs%variable_buoyforce .or. cs%first_call_set_forcing) .and. &
350  (.not.cs%adiabatic)) then
351  call set_net_mass_forcing(fluxes, forces, g, us)
352  endif
353 
354  cs%first_call_set_forcing = .false.
355 
356  call cpu_clock_end(id_clock_forcing)
357  call calltree_leave("set_forcing")
358 

References mom_forcing_type::allocate_mech_forcing(), buoyancy_forcing_const(), buoyancy_forcing_from_data_override(), buoyancy_forcing_from_files(), buoyancy_forcing_linear(), buoyancy_forcing_zero(), mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), dumbbell_surface_forcing::dumbbell_buoyancy_forcing(), id_clock_forcing, wind_forcing_1gyre(), wind_forcing_2gyre(), wind_forcing_by_data_override(), wind_forcing_const(), wind_forcing_from_file(), and wind_forcing_gyres().

Referenced by mom_main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ surface_forcing_end()

subroutine mom_surface_forcing::surface_forcing_end ( type(surface_forcing_cs), pointer  CS,
type(forcing), intent(inout), optional  fluxes 
)
private

Deallocate memory associated with the surface forcing module.

Parameters
cspointer to control struct returned by a previous surface_forcing_init call
[in,out]fluxesA structure containing thermodynamic forcing fields

Definition at line 1821 of file MOM_surface_forcing.F90.

1821  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1822  !! a previous surface_forcing_init call
1823  type(forcing), optional, intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
1824 ! Arguments: CS - A pointer to the control structure returned by a previous
1825 ! call to surface_forcing_init, it will be deallocated here.
1826 ! (inout) fluxes - A structure containing pointers to any possible
1827 ! forcing fields. Unused fields have NULL ptrs.
1828 
1829  if (present(fluxes)) call deallocate_forcing_type(fluxes)
1830 
1831 !#CTRL# call controlled_forcing_end(CS%ctrl_forcing_CSp)
1832 
1833  if (associated(cs)) deallocate(cs)
1834  cs => null()
1835 

◆ surface_forcing_init()

subroutine, public mom_surface_forcing::surface_forcing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(surface_forcing_cs), pointer  CS,
type(tracer_flow_control_cs), pointer  tracer_flow_CSp 
)

Initialize the surface forcing module.

Parameters
[in]timeThe current model time
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters
[in,out]diagstructure used to regulate diagnostic output
cspointer to control struct returned by a previous surface_forcing_init call
tracer_flow_cspForcing for tracers?

Definition at line 1398 of file MOM_surface_forcing.F90.

1398  type(time_type), intent(in) :: Time !< The current model time
1399  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1400  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1401  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1402  type(diag_ctrl), target, intent(inout) :: diag !< structure used to regulate diagnostic output
1403  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
1404  !! a previous surface_forcing_init call
1405  type(tracer_flow_control_CS), pointer :: tracer_flow_CSp !< Forcing for tracers?
1406 
1407  ! Local variables
1408  type(directories) :: dirs
1409  logical :: new_sim
1410  type(time_type) :: Time_frc
1411  ! This include declares and sets the variable "version".
1412 # include "version_variable.h"
1413  real :: flux_const_default ! The unscaled value of FLUXCONST [m day-1]
1414  logical :: default_2018_answers
1415  character(len=40) :: mdl = "MOM_surface_forcing" ! This module's name.
1416  character(len=200) :: filename, gust_file ! The name of the gustiness input file.
1417 
1418  if (associated(cs)) then
1419  call mom_error(warning, "surface_forcing_init called with an associated "// &
1420  "control structure.")
1421  return
1422  endif
1423  allocate(cs)
1424 
1425  id_clock_forcing=cpu_clock_id('(Ocean surface forcing)', grain=clock_module)
1426  call cpu_clock_begin(id_clock_forcing)
1427 
1428  cs%diag => diag
1429  if (associated(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
1430 
1431  ! Read all relevant parameters and write them to the model log.
1432  call log_version(param_file, mdl, version, '')
1433  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
1434  "If true, Temperature and salinity are used as state "//&
1435  "variables.", default=.true.)
1436  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, &
1437  "The directory in which all input files are found.", &
1438  default=".")
1439  cs%inputdir = slasher(cs%inputdir)
1440 
1441  call get_param(param_file, mdl, "ADIABATIC", cs%adiabatic, &
1442  "There are no diapycnal mass fluxes if ADIABATIC is "//&
1443  "true. This assumes that KD = KDML = 0.0 and that "//&
1444  "there is no buoyancy forcing, but makes the model "//&
1445  "faster by eliminating subroutine calls.", default=.false.)
1446  call get_param(param_file, mdl, "VARIABLE_WINDS", cs%variable_winds, &
1447  "If true, the winds vary in time after the initialization.", &
1448  default=.true.)
1449  call get_param(param_file, mdl, "VARIABLE_BUOYFORCE", cs%variable_buoyforce, &
1450  "If true, the buoyancy forcing varies in time after the "//&
1451  "initialization of the model.", default=.true.)
1452 
1453  call get_param(param_file, mdl, "BUOY_CONFIG", cs%buoy_config, &
1454  "The character string that indicates how buoyancy forcing "//&
1455  "is specified. Valid options include (file), (zero), "//&
1456  "(linear), (USER), (BFB) and (NONE).", fail_if_missing=.true.)
1457  if (trim(cs%buoy_config) == "file") then
1458  call get_param(param_file, mdl, "ARCHAIC_OMIP_FORCING_FILE", cs%archaic_OMIP_file, &
1459  "If true, use the forcing variable decomposition from "//&
1460  "the old German OMIP prescription that predated CORE. If "//&
1461  "false, use the variable groupings available from MOM "//&
1462  "output diagnostics of forcing variables.", default=.true.)
1463  if (cs%archaic_OMIP_file) then
1464  call get_param(param_file, mdl, "LONGWAVEDOWN_FILE", cs%longwave_file, &
1465  "The file with the downward longwave heat flux, in "//&
1466  "variable lwdn_sfc.", fail_if_missing=.true.)
1467  call get_param(param_file, mdl, "LONGWAVEUP_FILE", cs%longwaveup_file, &
1468  "The file with the upward longwave heat flux, in "//&
1469  "variable lwup_sfc.", fail_if_missing=.true.)
1470  call get_param(param_file, mdl, "EVAPORATION_FILE", cs%evaporation_file, &
1471  "The file with the evaporative moisture flux, in "//&
1472  "variable evap.", fail_if_missing=.true.)
1473  call get_param(param_file, mdl, "SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1474  "The file with the sensible heat flux, in "//&
1475  "variable shflx.", fail_if_missing=.true.)
1476  call get_param(param_file, mdl, "SHORTWAVEUP_FILE", cs%shortwaveup_file, &
1477  "The file with the upward shortwave heat flux.", &
1478  fail_if_missing=.true.)
1479  call get_param(param_file, mdl, "SHORTWAVEDOWN_FILE", cs%shortwave_file, &
1480  "The file with the downward shortwave heat flux.", &
1481  fail_if_missing=.true.)
1482  call get_param(param_file, mdl, "SNOW_FILE", cs%snow_file, &
1483  "The file with the downward frozen precip flux, in "//&
1484  "variable snow.", fail_if_missing=.true.)
1485  call get_param(param_file, mdl, "PRECIP_FILE", cs%rain_file, &
1486  "The file with the downward total precip flux, in "//&
1487  "variable precip.", fail_if_missing=.true.)
1488  call get_param(param_file, mdl, "FRESHDISCHARGE_FILE", cs%runoff_file, &
1489  "The file with the fresh and frozen runoff/calving fluxes, "//&
1490  "invariables disch_w and disch_s.", fail_if_missing=.true.)
1491 
1492  ! These variable names are hard-coded, per the archaic OMIP conventions.
1493  cs%latentheat_file = cs%evaporation_file ; cs%latent_var = "evap"
1494  cs%LW_var = "lwdn_sfc"; cs%SW_var = "swdn_sfc"; cs%sens_var = "shflx"
1495  cs%evap_var = "evap"; cs%rain_var = "precip"; cs%snow_var = "snow"
1496  cs%lrunoff_var = "disch_w"; cs%frunoff_var = "disch_s"
1497 
1498  else
1499  call get_param(param_file, mdl, "LONGWAVE_FILE", cs%longwave_file, &
1500  "The file with the longwave heat flux, in the variable "//&
1501  "given by LONGWAVE_FORCING_VAR.", fail_if_missing=.true.)
1502  call get_param(param_file, mdl, "LONGWAVE_FORCING_VAR", cs%LW_var, &
1503  "The variable with the longwave forcing field.", default="LW")
1504 
1505  call get_param(param_file, mdl, "SHORTWAVE_FILE", cs%shortwave_file, &
1506  "The file with the shortwave heat flux, in the variable "//&
1507  "given by SHORTWAVE_FORCING_VAR.", fail_if_missing=.true.)
1508  call get_param(param_file, mdl, "SHORTWAVE_FORCING_VAR", cs%SW_var, &
1509  "The variable with the shortwave forcing field.", default="SW")
1510 
1511  call get_param(param_file, mdl, "EVAPORATION_FILE", cs%evaporation_file, &
1512  "The file with the evaporative moisture flux, in the "//&
1513  "variable given by EVAP_FORCING_VAR.", fail_if_missing=.true.)
1514  call get_param(param_file, mdl, "EVAP_FORCING_VAR", cs%evap_var, &
1515  "The variable with the evaporative moisture flux.", &
1516  default="evap")
1517 
1518  call get_param(param_file, mdl, "LATENTHEAT_FILE", cs%latentheat_file, &
1519  "The file with the latent heat flux, in the variable "//&
1520  "given by LATENT_FORCING_VAR.", fail_if_missing=.true.)
1521  call get_param(param_file, mdl, "LATENT_FORCING_VAR", cs%latent_var, &
1522  "The variable with the latent heat flux.", default="latent")
1523 
1524  call get_param(param_file, mdl, "SENSIBLEHEAT_FILE", cs%sensibleheat_file, &
1525  "The file with the sensible heat flux, in the variable "//&
1526  "given by SENSIBLE_FORCING_VAR.", fail_if_missing=.true.)
1527  call get_param(param_file, mdl, "SENSIBLE_FORCING_VAR", cs%sens_var, &
1528  "The variable with the sensible heat flux.", default="sensible")
1529 
1530  call get_param(param_file, mdl, "RAIN_FILE", cs%rain_file, &
1531  "The file with the liquid precipitation flux, in the "//&
1532  "variable given by RAIN_FORCING_VAR.", fail_if_missing=.true.)
1533  call get_param(param_file, mdl, "RAIN_FORCING_VAR", cs%rain_var, &
1534  "The variable with the liquid precipitation flux.", &
1535  default="liq_precip")
1536  call get_param(param_file, mdl, "SNOW_FILE", cs%snow_file, &
1537  "The file with the frozen precipitation flux, in the "//&
1538  "variable given by SNOW_FORCING_VAR.", fail_if_missing=.true.)
1539  call get_param(param_file, mdl, "SNOW_FORCING_VAR", cs%snow_var, &
1540  "The variable with the frozen precipitation flux.", &
1541  default="froz_precip")
1542 
1543  call get_param(param_file, mdl, "RUNOFF_FILE", cs%runoff_file, &
1544  "The file with the fresh and frozen runoff/calving "//&
1545  "fluxes, in variables given by LIQ_RUNOFF_FORCING_VAR "//&
1546  "and FROZ_RUNOFF_FORCING_VAR.", fail_if_missing=.true.)
1547  call get_param(param_file, mdl, "LIQ_RUNOFF_FORCING_VAR", cs%lrunoff_var, &
1548  "The variable with the liquid runoff flux.", &
1549  default="liq_runoff")
1550  call get_param(param_file, mdl, "FROZ_RUNOFF_FORCING_VAR", cs%frunoff_var, &
1551  "The variable with the frozen runoff flux.", &
1552  default="froz_runoff")
1553  endif
1554 
1555  call get_param(param_file, mdl, "SSTRESTORE_FILE", cs%SSTrestore_file, &
1556  "The file with the SST toward which to restore in the "//&
1557  "variable given by SST_RESTORE_VAR.", fail_if_missing=.true.)
1558  call get_param(param_file, mdl, "SALINITYRESTORE_FILE", cs%salinityrestore_file, &
1559  "The file with the surface salinity toward which to "//&
1560  "restore in the variable given by SSS_RESTORE_VAR.", &
1561  fail_if_missing=.true.)
1562 
1563  if (cs%archaic_OMIP_file) then
1564  cs%SST_restore_var = "TEMP" ; cs%SSS_restore_var = "SALT"
1565  else
1566  call get_param(param_file, mdl, "SST_RESTORE_VAR", cs%SST_restore_var, &
1567  "The variable with the SST toward which to restore.", &
1568  default="SST")
1569  call get_param(param_file, mdl, "SSS_RESTORE_VAR", cs%SSS_restore_var, &
1570  "The variable with the SSS toward which to restore.", &
1571  default="SSS")
1572  endif
1573 
1574  ! Add inputdir to the file names.
1575  cs%shortwave_file = trim(cs%inputdir)//trim(cs%shortwave_file)
1576  cs%longwave_file = trim(cs%inputdir)//trim(cs%longwave_file)
1577  cs%sensibleheat_file = trim(cs%inputdir)//trim(cs%sensibleheat_file)
1578  cs%latentheat_file = trim(cs%inputdir)//trim(cs%latentheat_file)
1579  cs%evaporation_file = trim(cs%inputdir)//trim(cs%evaporation_file)
1580  cs%snow_file = trim(cs%inputdir)//trim(cs%snow_file)
1581  cs%rain_file = trim(cs%inputdir)//trim(cs%rain_file)
1582  cs%runoff_file = trim(cs%inputdir)//trim(cs%runoff_file)
1583 
1584  cs%shortwaveup_file = trim(cs%inputdir)//trim(cs%shortwaveup_file)
1585  cs%longwaveup_file = trim(cs%inputdir)//trim(cs%longwaveup_file)
1586 
1587  cs%SSTrestore_file = trim(cs%inputdir)//trim(cs%SSTrestore_file)
1588  cs%salinityrestore_file = trim(cs%inputdir)//trim(cs%salinityrestore_file)
1589  elseif (trim(cs%buoy_config) == "const") then
1590  call get_param(param_file, mdl, "SENSIBLE_HEAT_FLUX", cs%constantHeatForcing, &
1591  "A constant heat forcing (positive into ocean) applied "//&
1592  "through the sensible heat flux field. ", &
1593  units='W/m2', fail_if_missing=.true.)
1594  endif
1595  call get_param(param_file, mdl, "WIND_CONFIG", cs%wind_config, &
1596  "The character string that indicates how wind forcing "//&
1597  "is specified. Valid options include (file), (2gyre), "//&
1598  "(1gyre), (gyres), (zero), and (USER).", fail_if_missing=.true.)
1599  if (trim(cs%wind_config) == "file") then
1600  call get_param(param_file, mdl, "WIND_FILE", cs%wind_file, &
1601  "The file in which the wind stresses are found in "//&
1602  "variables STRESS_X and STRESS_Y.", fail_if_missing=.true.)
1603  call get_param(param_file, mdl, "WINDSTRESS_X_VAR",cs%stress_x_var, &
1604  "The name of the x-wind stress variable in WIND_FILE.", &
1605  default="STRESS_X")
1606  call get_param(param_file, mdl, "WINDSTRESS_Y_VAR", cs%stress_y_var, &
1607  "The name of the y-wind stress variable in WIND_FILE.", &
1608  default="STRESS_Y")
1609  call get_param(param_file, mdl, "WINDSTRESS_STAGGER",cs%wind_stagger, &
1610  "A character indicating how the wind stress components "//&
1611  "are staggered in WIND_FILE. This may be A or C for now.", &
1612  default="A")
1613  call get_param(param_file, mdl, "WINDSTRESS_SCALE", cs%wind_scale, &
1614  "A value by which the wind stresses in WIND_FILE are rescaled.", &
1615  default=1.0, units="nondim")
1616  call get_param(param_file, mdl, "USTAR_FORCING_VAR", cs%ustar_var, &
1617  "The name of the friction velocity variable in WIND_FILE "//&
1618  "or blank to get ustar from the wind stresses plus the "//&
1619  "gustiness.", default=" ", units="nondim")
1620  cs%wind_file = trim(cs%inputdir) // trim(cs%wind_file)
1621  endif
1622  if (trim(cs%wind_config) == "gyres") then
1623  call get_param(param_file, mdl, "TAUX_CONST", cs%gyres_taux_const, &
1624  "With the gyres wind_config, the constant offset in the "//&
1625  "zonal wind stress profile: "//&
1626  " A in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1627  units="Pa", default=0.0)
1628  call get_param(param_file, mdl, "TAUX_SIN_AMP",cs%gyres_taux_sin_amp, &
1629  "With the gyres wind_config, the sine amplitude in the "//&
1630  "zonal wind stress profile: "//&
1631  " B in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1632  units="Pa", default=0.0)
1633  call get_param(param_file, mdl, "TAUX_COS_AMP",cs%gyres_taux_cos_amp, &
1634  "With the gyres wind_config, the cosine amplitude in "//&
1635  "the zonal wind stress profile: "//&
1636  " C in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1637  units="Pa", default=0.0)
1638  call get_param(param_file, mdl, "TAUX_N_PIS",cs%gyres_taux_n_pis, &
1639  "With the gyres wind_config, the number of gyres in "//&
1640  "the zonal wind stress profile: "//&
1641  " n in taux = A + B*sin(n*pi*y/L) + C*cos(n*pi*y/L).", &
1642  units="nondim", default=0.0)
1643  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1644  "This sets the default value for the various _2018_ANSWERS parameters.", &
1645  default=.true.)
1646  call get_param(param_file, mdl, "WIND_GYRES_2018_ANSWERS", cs%answers_2018, &
1647  "If true, use the order of arithmetic and expressions that recover the answers "//&
1648  "from the end of 2018. Otherwise, use expressions for the gyre friction velocities "//&
1649  "that are rotationally invariant and more likely to be the same between compilers.", &
1650  default=default_2018_answers)
1651  else
1652  cs%answers_2018 = .false.
1653  endif
1654  if ((trim(cs%wind_config) == "2gyre") .or. &
1655  (trim(cs%wind_config) == "1gyre") .or. &
1656  (trim(cs%wind_config) == "gyres") .or. &
1657  (trim(cs%buoy_config) == "linear")) then
1658  cs%south_lat = g%south_lat
1659  cs%len_lat = g%len_lat
1660  endif
1661  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
1662  "The mean ocean density used with BOUSSINESQ true to "//&
1663  "calculate accelerations and the mass for conservation "//&
1664  "properties, or with BOUSSINSEQ false to convert some "//&
1665  "parameters from vertical units of m to kg m-2.", &
1666  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
1667  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
1668  "If true, the buoyancy fluxes drive the model back "//&
1669  "toward some specified surface state with a rate "//&
1670  "given by FLUXCONST.", default= .false.)
1671  call get_param(param_file, mdl, "LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
1672  "The latent heat of fusion.", default=hlf, &
1673  units="J/kg", scale=us%R_to_kg_m3*us%Z_to_m*us%s_to_T)
1674  call get_param(param_file, mdl, "LATENT_HEAT_VAPORIZATION", cs%latent_heat_vapor, &
1675  "The latent heat of fusion.", units="J/kg", default=hlv)
1676  if (cs%restorebuoy) then
1677  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
1678  "The constant that relates the restoring surface fluxes "//&
1679  "to the relative surface anomalies (akin to a piston "//&
1680  "velocity). Note the non-MKS units.", &
1681  units="m day-1", scale=us%m_to_Z*us%T_to_s, &
1682  fail_if_missing=.true., unscaled=flux_const_default)
1683 
1684  if (cs%use_temperature) then
1685  call get_param(param_file, mdl, "FLUXCONST_T", cs%Flux_const_T, &
1686  "The constant that relates the restoring surface temperature "//&
1687  "flux to the relative surface anomaly (akin to a piston "//&
1688  "velocity). Note the non-MKS units.", &
1689  units="m day-1", scale=1.0, & ! scale=US%m_to_Z*US%T_to_s,
1690  default=flux_const_default)
1691  call get_param(param_file, mdl, "FLUXCONST_S", cs%Flux_const_S, &
1692  "The constant that relates the restoring surface salinity "//&
1693  "flux to the relative surface anomaly (akin to a piston "//&
1694  "velocity). Note the non-MKS units.", &
1695  units="m day-1", scale=us%m_to_Z*us%T_to_s, &
1696  default=flux_const_default)
1697  endif
1698 
1699  !### Convert flux constants from m day-1 to m s-1. Folding these into the scaling
1700  ! factors above could change a division into a multiply by a reciprocal, which could
1701  ! change answers at the level of roundoff.
1702  cs%Flux_const = cs%Flux_const / 86400.0
1703  cs%Flux_const_T = cs%Flux_const_T / 86400.0
1704  cs%Flux_const_S = cs%Flux_const_S / 86400.0
1705 
1706  if (trim(cs%buoy_config) == "linear") then
1707  call get_param(param_file, mdl, "SST_NORTH", cs%T_north, &
1708  "With buoy_config linear, the sea surface temperature "//&
1709  "at the northern end of the domain toward which to "//&
1710  "to restore.", units="deg C", default=0.0)
1711  call get_param(param_file, mdl, "SST_SOUTH", cs%T_south, &
1712  "With buoy_config linear, the sea surface temperature "//&
1713  "at the southern end of the domain toward which to "//&
1714  "to restore.", units="deg C", default=0.0)
1715  call get_param(param_file, mdl, "SSS_NORTH", cs%S_north, &
1716  "With buoy_config linear, the sea surface salinity "//&
1717  "at the northern end of the domain toward which to "//&
1718  "to restore.", units="PSU", default=35.0)
1719  call get_param(param_file, mdl, "SSS_SOUTH", cs%S_south, &
1720  "With buoy_config linear, the sea surface salinity "//&
1721  "at the southern end of the domain toward which to "//&
1722  "to restore.", units="PSU", default=35.0)
1723  endif
1724  endif
1725  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
1726  "The gravitational acceleration of the Earth.", &
1727  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
1728 
1729  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
1730  "The background gustiness in the winds.", &
1731  units="Pa", default=0.02, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
1732  call get_param(param_file, mdl, "READ_GUST_2D", cs%read_gust_2d, &
1733  "If true, use a 2-dimensional gustiness supplied from "//&
1734  "an input file", default=.false.)
1735  if (cs%read_gust_2d) then
1736  call get_param(param_file, mdl, "GUST_2D_FILE", gust_file, &
1737  "The file in which the wind gustiness is found in "//&
1738  "variable gustiness.", fail_if_missing=.true.)
1739  call safe_alloc_ptr(cs%gust,g%isd,g%ied,g%jsd,g%jed)
1740  filename = trim(cs%inputdir) // trim(gust_file)
1741  call mom_read_data(filename,'gustiness',cs%gust,g%domain, timelevel=1, &
1742  scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z) ! units in file should be Pa
1743  endif
1744 
1745 ! All parameter settings are now known.
1746 
1747  if (trim(cs%wind_config) == "USER" .or. trim(cs%buoy_config) == "USER" ) then
1748  call user_surface_forcing_init(time, g, us, param_file, diag, cs%user_forcing_CSp)
1749  elseif (trim(cs%buoy_config) == "BFB" ) then
1750  call bfb_surface_forcing_init(time, g, us, param_file, diag, cs%BFB_forcing_CSp)
1751  elseif (trim(cs%buoy_config) == "dumbbell" ) then
1752  call dumbbell_surface_forcing_init(time, g, us, param_file, diag, cs%dumbbell_forcing_CSp)
1753  elseif (trim(cs%wind_config) == "MESO" .or. trim(cs%buoy_config) == "MESO" ) then
1754  call meso_surface_forcing_init(time, g, us, param_file, diag, cs%MESO_forcing_CSp)
1755  elseif (trim(cs%wind_config) == "Neverland") then
1756  call neverland_surface_forcing_init(time, g, us, param_file, diag, cs%Neverland_forcing_CSp)
1757  elseif (trim(cs%wind_config) == "ideal_hurr" .or.&
1758  trim(cs%wind_config) == "SCM_ideal_hurr") then
1759  call idealized_hurricane_wind_init(time, g, us, param_file, cs%idealized_hurricane_CSp)
1760  elseif (trim(cs%wind_config) == "const") then
1761  call get_param(param_file, mdl, "CONST_WIND_TAUX", cs%tau_x0, &
1762  "With wind_config const, this is the constant zonal "//&
1763  "wind-stress", units="Pa", fail_if_missing=.true.)
1764  call get_param(param_file, mdl, "CONST_WIND_TAUY", cs%tau_y0, &
1765  "With wind_config const, this is the constant meridional "//&
1766  "wind-stress", units="Pa", fail_if_missing=.true.)
1767  elseif (trim(cs%wind_config) == "SCM_CVmix_tests" .or. &
1768  trim(cs%buoy_config) == "SCM_CVmix_tests") then
1769  call scm_cvmix_tests_surface_forcing_init(time, g, param_file, cs%SCM_CVmix_tests_CSp)
1770  cs%SCM_CVmix_tests_CSp%Rho0 = us%R_to_kg_m3*cs%Rho0 !copy reference density for pass
1771  endif
1772 
1773  call register_forcing_type_diags(time, diag, us, cs%use_temperature, cs%handles)
1774 
1775  ! Set up any restart fields associated with the forcing.
1776  call restart_init(param_file, cs%restart_CSp, "MOM_forcing.res")
1777 !#CTRL# call register_ctrl_forcing_restarts(G, param_file, CS%ctrl_forcing_CSp, &
1778 !#CTRL# CS%restart_CSp)
1779  call restart_init_end(cs%restart_CSp)
1780 
1781  if (associated(cs%restart_CSp)) then
1782  call get_mom_input(dirs=dirs)
1783 
1784  new_sim = .false.
1785  if ((dirs%input_filename(1:1) == 'n') .and. &
1786  (len_trim(dirs%input_filename) == 1)) new_sim = .true.
1787  if (.not.new_sim) then
1788  call restore_state(dirs%input_filename, dirs%restart_input_dir, time_frc, &
1789  g, cs%restart_CSp)
1790  endif
1791  endif
1792 
1793  ! Determine how many time levels are in each forcing variable.
1794  if (trim(cs%buoy_config) == "file") then
1795  cs%SW_nlev = num_timelevels(cs%shortwave_file, cs%SW_var, min_dims=3)
1796  cs%LW_nlev = num_timelevels(cs%longwave_file, cs%LW_var, min_dims=3)
1797  cs%latent_nlev = num_timelevels(cs%latentheat_file, cs%latent_var, 3)
1798  cs%sens_nlev = num_timelevels(cs%sensibleheat_file, cs%sens_var, min_dims=3)
1799 
1800  cs%evap_nlev = num_timelevels(cs%evaporation_file, cs%evap_var, min_dims=3)
1801  cs%precip_nlev = num_timelevels(cs%rain_file, cs%rain_var, min_dims=3)
1802  cs%runoff_nlev = num_timelevels(cs%runoff_file, cs%lrunoff_var, 3)
1803 
1804  cs%SST_nlev = num_timelevels(cs%SSTrestore_file, cs%SST_restore_var, 3)
1805  cs%SSS_nlev = num_timelevels(cs%salinityrestore_file, cs%SSS_restore_var, 3)
1806  endif
1807 
1808  if (trim(cs%wind_config) == "file") &
1809  cs%wind_nlev = num_timelevels(cs%wind_file, cs%stress_x_var, min_dims=3)
1810 
1811 !#CTRL# call controlled_forcing_init(Time, G, param_file, diag, CS%ctrl_forcing_CSp)
1812 
1813  call user_revise_forcing_init(param_file, cs%urf_CS)
1814 
1815  call cpu_clock_end(id_clock_forcing)

References bfb_surface_forcing::bfb_surface_forcing_init(), mom_get_input::get_mom_input(), id_clock_forcing, meso_surface_forcing::meso_surface_forcing_init(), neverland_surface_forcing::neverland_surface_forcing_init(), mom_io::num_timelevels(), mom_restart::restart_init_end(), mom_restart::restore_state(), and user_surface_forcing::user_surface_forcing_init().

Referenced by mom_main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_1gyre()

subroutine mom_surface_forcing::wind_forcing_1gyre ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Sets the surface wind stresses to set up a single idealized gyre.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]dayThe time of the fluxes
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 445 of file MOM_surface_forcing.F90.

445  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
446  !! describe the surface state of the ocean.
447  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
448  type(time_type), intent(in) :: day !< The time of the fluxes
449  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
450  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
451  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
452  !! a previous surface_forcing_init call
453  ! Local variables
454  real :: PI
455  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
456 
457  call calltree_enter("wind_forcing_1gyre, MOM_surface_forcing.F90")
458  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
459  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
460 
461  ! set the steady surface wind stresses, in units of Pa.
462  pi = 4.0*atan(1.0)
463 
464  do j=js,je ; do i=is-1,ieq
465  forces%taux(i,j) = -0.2*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
466  cos(pi*(g%geoLatCu(i,j)-cs%South_lat)/cs%len_lat)
467  enddo ; enddo
468 
469  do j=js-1,jeq ; do i=is,ie
470  forces%tauy(i,j) = 0.0
471  enddo ; enddo
472 
473  call calltree_leave("wind_forcing_1gyre")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_2gyre()

subroutine mom_surface_forcing::wind_forcing_2gyre ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Sets the surface wind stresses to set up two idealized gyres.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]dayThe time of the fluxes
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 411 of file MOM_surface_forcing.F90.

411  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
412  !! describe the surface state of the ocean.
413  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
414  type(time_type), intent(in) :: day !< The time of the fluxes
415  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
416  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
417  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
418  !! a previous surface_forcing_init call
419  ! Local variables
420  real :: PI
421  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
422 
423  call calltree_enter("wind_forcing_2gyre, MOM_surface_forcing.F90")
424  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
425  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
426 
427  !set the steady surface wind stresses, in units of Pa.
428  pi = 4.0*atan(1.0)
429 
430  do j=js,je ; do i=is-1,ieq
431  forces%taux(i,j) = 0.1*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
432  (1.0 - cos(2.0*pi*(g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat))
433  enddo ; enddo
434 
435  do j=js-1,jeq ; do i=is,ie
436  forces%tauy(i,j) = 0.0
437  enddo ; enddo
438 
439  call calltree_leave("wind_forcing_2gyre")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_by_data_override()

subroutine mom_surface_forcing::wind_forcing_by_data_override ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private
Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]dayThe time of the fluxes
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 686 of file MOM_surface_forcing.F90.

686  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
687  !! describe the surface state of the ocean.
688  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
689  type(time_type), intent(in) :: day !< The time of the fluxes
690  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
691  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
692  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
693  !! a previous surface_forcing_init call
694  ! Local variables
695  real :: temp_x(SZI_(G),SZJ_(G)) ! Pseudo-zonal and psuedo-meridional
696  real :: temp_y(SZI_(G),SZJ_(G)) ! wind stresses at h-points [Pa].
697  real :: temp_ustar(SZI_(G),SZJ_(G)) ! ustar [m s-1] (not rescaled).
698  real :: Pa_conversion ! A unit conversion factor from Pa to the internal units [R Z L T-2 Pa-1 ~> 1]
699  integer :: i, j, is_in, ie_in, js_in, je_in
700  logical :: read_uStar
701 
702  call calltree_enter("wind_forcing_by_data_override, MOM_surface_forcing.F90")
703 
704  if (.not.cs%dataOverrideIsInitialized) then
705  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true., press=.true.)
706  call data_override_init(ocean_domain_in=g%Domain%mpp_domain)
707  cs%dataOverrideIsInitialized = .true.
708  endif
709 
710  is_in = g%isc - g%isd + 1 ; ie_in = g%iec - g%isd + 1
711  js_in = g%jsc - g%jsd + 1 ; je_in = g%jec - g%jsd + 1
712  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
713 
714  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
715  call data_override('OCN', 'taux', temp_x, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
716  call data_override('OCN', 'tauy', temp_y, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
717  call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
718  ! Ignore CS%wind_scale when using data_override ?????
719  do j=g%jsc,g%jec ; do i=g%isc-1,g%IecB
720  forces%taux(i,j) = pa_conversion * 0.5 * (temp_x(i,j) + temp_x(i+1,j))
721  enddo ; enddo
722  do j=g%jsc-1,g%JecB ; do i=g%isc,g%iec
723  forces%tauy(i,j) = pa_conversion * 0.5 * (temp_y(i,j) + temp_y(i,j+1))
724  enddo ; enddo
725 
726  read_ustar = (len_trim(cs%ustar_var) > 0) ! Need better control higher up ????
727  if (read_ustar) then
728  do j=g%jsc,g%jec ; do i=g%isc,g%iec ; temp_ustar(i,j) = us%Z_to_m*us%s_to_T*forces%ustar(i,j) ; enddo ; enddo
729  call data_override('OCN', 'ustar', temp_ustar, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
730  do j=g%jsc,g%jec ; do i=g%isc,g%iec ; forces%ustar(i,j) = us%m_to_Z*us%T_to_s*temp_ustar(i,j) ; enddo ; enddo
731  else
732  if (cs%read_gust_2d) then
733  call data_override('OCN', 'gust', cs%gust, day, is_in=is_in, ie_in=ie_in, js_in=js_in, je_in=je_in)
734  do j=g%jsc,g%jec ; do i=g%isc,g%iec
735  forces%ustar(i,j) = sqrt((pa_conversion * sqrt(temp_x(i,j)*temp_x(i,j) + &
736  temp_y(i,j)*temp_y(i,j)) + cs%gust(i,j)) * us%L_to_Z / cs%Rho0)
737  enddo ; enddo
738  else
739  do j=g%jsc,g%jec ; do i=g%isc,g%iec
740  forces%ustar(i,j) = sqrt(us%L_to_Z * (pa_conversion*sqrt(temp_x(i,j)*temp_x(i,j) + &
741  temp_y(i,j)*temp_y(i,j))/cs%Rho0 + cs%gust_const/cs%Rho0 ))
742  enddo ; enddo
743  endif
744  endif
745 
746  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
747 ! call pass_var(forces%ustar, G%Domain, To_All) Not needed ?????
748 
749  call calltree_leave("wind_forcing_by_data_override")

References mom_forcing_type::allocate_mech_forcing(), mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_const()

subroutine mom_surface_forcing::wind_forcing_const ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
real, intent(in)  tau_x0,
real, intent(in)  tau_y0,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Sets the surface wind stresses to constant values.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]tau_x0The zonal wind stress [Pa]
[in]tau_y0The meridional wind stress [Pa]
[in]dayThe time of the fluxes
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 363 of file MOM_surface_forcing.F90.

363  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
364  !! describe the surface state of the ocean.
365  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
366  real, intent(in) :: tau_x0 !< The zonal wind stress [Pa]
367  real, intent(in) :: tau_y0 !< The meridional wind stress [Pa]
368  type(time_type), intent(in) :: day !< The time of the fluxes
369  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
370  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
371  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
372  !! a previous surface_forcing_init call
373  ! Local variables
374  real :: Pa_conversion ! A unit conversion factor from Pa to the internal units [R Z L T-2 Pa-1 ~> 1]
375  real :: mag_tau
376  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
377 
378  call calltree_enter("wind_forcing_const, MOM_surface_forcing.F90")
379  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
380  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
381  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
382 
383  !set steady surface wind stresses, in units of Pa.
384  !### mag_tau = US%kg_m3_to_R*US%m_s_to_L_T**2*US%L_to_Z * sqrt( tau_x0**2 + tau_y0**2)
385  mag_tau = pa_conversion * sqrt( tau_x0**2 + tau_y0**2)
386 
387  do j=js,je ; do i=is-1,ieq
388  forces%taux(i,j) = tau_x0 * pa_conversion
389  enddo ; enddo
390 
391  do j=js-1,jeq ; do i=is,ie
392  forces%tauy(i,j) = tau_y0 * pa_conversion
393  enddo ; enddo
394 
395  if (cs%read_gust_2d) then
396  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
397  forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust(i,j) ) / cs%Rho0 )
398  enddo ; enddo ; endif
399  else
400  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
401  forces%ustar(i,j) = sqrt( us%L_to_Z * ( mag_tau + cs%gust_const ) / cs%Rho0 )
402  enddo ; enddo ; endif
403  endif
404 
405  call calltree_leave("wind_forcing_const")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_from_file()

subroutine mom_surface_forcing::wind_forcing_from_file ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private
Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]dayThe time of the fluxes
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 531 of file MOM_surface_forcing.F90.

531  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
532  !! describe the surface state of the ocean.
533  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
534  type(time_type), intent(in) :: day !< The time of the fluxes
535  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
536  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
537  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
538  !! a previous surface_forcing_init call
539  ! Local variables
540  character(len=200) :: filename ! The name of the input file.
541  real :: temp_x(SZI_(G),SZJ_(G)) ! Pseudo-zonal and psuedo-meridional
542  real :: temp_y(SZI_(G),SZJ_(G)) ! wind stresses at h-points [R L Z T-1 ~> Pa].
543  real :: Pa_conversion ! A unit conversion factor from Pa to the internal wind stress
544  ! units [R Z L T-2 Pa-1 ~> 1]
545  integer :: time_lev_daily ! The time levels to read for fields with
546  integer :: time_lev_monthly ! daily and montly cycles.
547  integer :: time_lev ! The time level that is used for a field.
548  integer :: days, seconds
549  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
550  logical :: read_Ustar
551 
552  call calltree_enter("wind_forcing_from_file, MOM_surface_forcing.F90")
553  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
554  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
555  pa_conversion = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
556 
557  call get_time(day, seconds, days)
558  time_lev_daily = days - 365*floor(real(days) / 365.0)
559 
560  if (time_lev_daily < 31) then ; time_lev_monthly = 0
561  elseif (time_lev_daily < 59) then ; time_lev_monthly = 1
562  elseif (time_lev_daily < 90) then ; time_lev_monthly = 2
563  elseif (time_lev_daily < 120) then ; time_lev_monthly = 3
564  elseif (time_lev_daily < 151) then ; time_lev_monthly = 4
565  elseif (time_lev_daily < 181) then ; time_lev_monthly = 5
566  elseif (time_lev_daily < 212) then ; time_lev_monthly = 6
567  elseif (time_lev_daily < 243) then ; time_lev_monthly = 7
568  elseif (time_lev_daily < 273) then ; time_lev_monthly = 8
569  elseif (time_lev_daily < 304) then ; time_lev_monthly = 9
570  elseif (time_lev_daily < 334) then ; time_lev_monthly = 10
571  else ; time_lev_monthly = 11
572  endif
573 
574  time_lev_daily = time_lev_daily+1
575  time_lev_monthly = time_lev_monthly+1
576 
577  select case (cs%wind_nlev)
578  case (12) ; time_lev = time_lev_monthly
579  case (365) ; time_lev = time_lev_daily
580  case default ; time_lev = 1
581  end select
582 
583  if (time_lev /= cs%wind_last_lev) then
584  filename = trim(cs%wind_file)
585  read_ustar = (len_trim(cs%ustar_var) > 0)
586 ! if (is_root_pe()) &
587 ! write(*,'("Wind_forcing Reading time level ",I," last was ",I,".")')&
588 ! time_lev-1,CS%wind_last_lev-1
589  select case ( uppercase(cs%wind_stagger(1:1)) )
590  case ("A")
591  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
592  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
593  temp_x(:,:), temp_y(:,:), g%Domain, stagger=agrid, &
594  timelevel=time_lev, scale=pa_conversion)
595 
596  call pass_vector(temp_x, temp_y, g%Domain, to_all, agrid)
597  do j=js,je ; do i=is-1,ieq
598  forces%taux(i,j) = 0.5 * cs%wind_scale * (temp_x(i,j) + temp_x(i+1,j))
599  enddo ; enddo
600  do j=js-1,jeq ; do i=is,ie
601  forces%tauy(i,j) = 0.5 * cs%wind_scale * (temp_y(i,j) + temp_y(i,j+1))
602  enddo ; enddo
603 
604  if (.not.read_ustar) then
605  if (cs%read_gust_2d) then
606  do j=js,je ; do i=is,ie
607  forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
608  sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j))) * us%L_to_Z / cs%Rho0)
609  enddo ; enddo
610  else
611  do j=js,je ; do i=is,ie
612  forces%ustar(i,j) = sqrt(us%L_to_Z * (cs%gust_const/cs%Rho0 + &
613  sqrt(temp_x(i,j)*temp_x(i,j) + temp_y(i,j)*temp_y(i,j)) / cs%Rho0) )
614  enddo ; enddo
615  endif
616  endif
617  case ("C")
618  if (g%symmetric) then
619  if (.not.associated(g%Domain_aux)) call mom_error(fatal, &
620  " wind_forcing_from_file with C-grid input and symmetric memory "//&
621  " called with a non-associated auxiliary domain in the grid type.")
622  ! Read the data as though symmetric memory were not being used, and
623  ! then translate it appropriately.
624  temp_x(:,:) = 0.0 ; temp_y(:,:) = 0.0
625  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
626  temp_x(:,:), temp_y(:,:), &
627  g%Domain_aux, stagger=cgrid_ne, timelevel=time_lev, &
628  scale=pa_conversion)
629  do j=js,je ; do i=is,ie
630  forces%taux(i,j) = cs%wind_scale * temp_x(i,j)
631  forces%tauy(i,j) = cs%wind_scale * temp_y(i,j)
632  enddo ; enddo
633  call fill_symmetric_edges(forces%taux, forces%tauy, g%Domain, stagger=cgrid_ne)
634  else
635  call mom_read_vector(filename, cs%stress_x_var, cs%stress_y_var, &
636  forces%taux(:,:), forces%tauy(:,:), &
637  g%Domain, stagger=cgrid_ne, timelevel=time_lev, &
638  scale=pa_conversion)
639 
640  if (cs%wind_scale /= 1.0) then
641  do j=js,je ; do i=isq,ieq
642  forces%taux(i,j) = cs%wind_scale * forces%taux(i,j)
643  enddo ; enddo
644  do j=jsq,jeq ; do i=is,ie
645  forces%tauy(i,j) = cs%wind_scale * forces%tauy(i,j)
646  enddo ; enddo
647  endif
648  endif
649 
650  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
651  if (.not.read_ustar) then
652  if (cs%read_gust_2d) then
653  do j=js, je ; do i=is, ie
654  forces%ustar(i,j) = sqrt((cs%gust(i,j) + &
655  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
656  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * us%L_to_Z / cs%Rho0 )
657  enddo ; enddo
658  else
659  do j=js, je ; do i=is, ie
660  forces%ustar(i,j) = sqrt(us%L_to_Z * ( (cs%gust_const/cs%Rho0) + &
661  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
662  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2)))/cs%Rho0))
663  enddo ; enddo
664  endif
665  endif
666  case default
667  call mom_error(fatal, "wind_forcing_from_file: Unrecognized stagger "//&
668  trim(cs%wind_stagger)//" is not 'A' or 'C'.")
669  end select
670 
671  if (read_ustar) then
672  call mom_read_data(filename, cs%Ustar_var, forces%ustar(:,:), &
673  g%Domain, timelevel=time_lev, scale=us%m_to_Z*us%T_to_s)
674  endif
675 
676  cs%wind_last_lev = time_lev
677 
678  endif ! time_lev /= CS%wind_last_lev
679 
680  call calltree_leave("wind_forcing_from_file")

References mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), and mom_string_functions::uppercase().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wind_forcing_gyres()

subroutine mom_surface_forcing::wind_forcing_gyres ( type(surface), intent(inout)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(surface_forcing_cs), pointer  CS 
)
private

Sets the surface wind stresses to set up idealized gyres.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]forcesA structure with the driving mechanical forces
[in]dayThe time of the fluxes
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
cspointer to control struct returned by a previous surface_forcing_init call

Definition at line 478 of file MOM_surface_forcing.F90.

478  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
479  !! describe the surface state of the ocean.
480  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
481  type(time_type), intent(in) :: day !< The time of the fluxes
482  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
483  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
484  type(surface_forcing_CS), pointer :: CS !< pointer to control struct returned by
485  !! a previous surface_forcing_init call
486  ! Local variables
487  real :: PI, y, I_rho
488  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
489 
490  call calltree_enter("wind_forcing_gyres, MOM_surface_forcing.F90")
491  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
492  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
493 
494  ! steady surface wind stresses [Pa]
495  pi = 4.0*atan(1.0)
496 
497  do j=js-1,je+1 ; do i=is-1,ieq
498  y = (g%geoLatCu(i,j)-cs%South_lat) / cs%len_lat
499  forces%taux(i,j) = us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z * &
500  (cs%gyres_taux_const + &
501  ( cs%gyres_taux_sin_amp*sin(cs%gyres_taux_n_pis*pi*y) &
502  + cs%gyres_taux_cos_amp*cos(cs%gyres_taux_n_pis*pi*y) ))
503  enddo ; enddo
504 
505  do j=js-1,jeq ; do i=is-1,ie+1
506  forces%tauy(i,j) = 0.0
507  enddo ; enddo
508 
509  ! set the friction velocity
510  if (cs%answers_2018) then
511  do j=js,je ; do i=is,ie
512  forces%ustar(i,j) = sqrt(us%L_to_Z * ((cs%gust_const/cs%Rho0) + &
513  sqrt(0.5*(forces%tauy(i,j-1)*forces%tauy(i,j-1) + forces%tauy(i,j)*forces%tauy(i,j) + &
514  forces%taux(i-1,j)*forces%taux(i-1,j) + forces%taux(i,j)*forces%taux(i,j)))/cs%Rho0) )
515  enddo ; enddo
516  else
517  i_rho = us%L_to_Z / cs%Rho0
518  do j=js,je ; do i=is,ie
519  forces%ustar(i,j) = sqrt( (cs%gust_const + &
520  sqrt(0.5*((forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2) + &
521  (forces%taux(i-1,j)**2 + forces%taux(i,j)**2))) ) * i_rho )
522  enddo ; enddo
523  endif
524 
525  call calltree_leave("wind_forcing_gyres")

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by set_forcing().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ id_clock_forcing

integer mom_surface_forcing::id_clock_forcing

A CPU time clock.

Definition at line 212 of file MOM_surface_forcing.F90.

212 integer :: id_clock_forcing !< A CPU time clock

Referenced by set_forcing(), and surface_forcing_init().