MOM6
user_surface_forcing Module Reference

Detailed Description

Template for user to code up surface forcing.

Data Types

type  user_surface_forcing_cs
 This control structure should be used to store any run-time variables associated with the user-specified forcing. More...
 

Functions/Subroutines

subroutine, public user_wind_forcing (sfc_state, forces, day, G, US, CS)
 This subroutine sets the surface wind stresses, forcestaux and forcestauy, in [R Z L T-2 ~> Pa]. These are the stresses in the direction of the model grid (i.e. the same direction as the u- and v- velocities). More...
 
subroutine, public user_buoyancy_forcing (sfc_state, fluxes, day, dt, G, US, CS)
 This subroutine specifies the current surface fluxes of buoyancy or temperature and fresh water. It may also be modified to add surface fluxes of user provided tracers. More...
 
subroutine, public user_surface_forcing_init (Time, G, US, param_file, diag, CS)
 This subroutine initializes the USER_surface_forcing module. More...
 

Function/Subroutine Documentation

◆ user_buoyancy_forcing()

subroutine, public user_surface_forcing::user_buoyancy_forcing ( 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(user_surface_forcing_cs), pointer  CS 
)

This subroutine specifies the current surface fluxes of buoyancy or temperature and fresh water. It may also be modified to add surface fluxes of user provided tracers.

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
csA pointer to the control structure returned by a previous call to user_surface_forcing_init

Definition at line 103 of file user_surface_forcing.F90.

103  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
104  !! describe the surface state of the ocean.
105  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
106  type(time_type), intent(in) :: day !< The time of the fluxes
107  real, intent(in) :: dt !< The amount of time over which
108  !! the fluxes apply [s]
109  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
110  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
111  type(user_surface_forcing_CS), pointer :: CS !< A pointer to the control structure returned
112  !! by a previous call to user_surface_forcing_init
113 
114 ! This subroutine specifies the current surface fluxes of buoyancy or
115 ! temperature and fresh water. It may also be modified to add
116 ! surface fluxes of user provided tracers.
117 
118 ! When temperature is used, there are long list of fluxes that need to be
119 ! set - essentially the same as for a full coupled model, but most of these
120 ! can be simply set to zero. The net fresh water flux should probably be
121 ! set in fluxes%evap and fluxes%lprec, with any salinity restoring
122 ! appearing in fluxes%vprec, and the other water flux components
123 ! (fprec, lrunoff and frunoff) left as arrays full of zeros.
124 ! Evap is usually negative and precip is usually positive. All heat fluxes
125 ! are in W m-2 and positive for heat going into the ocean. All fresh water
126 ! fluxes are in [R Z T-1 ~> kg m-2 s-1] and positive for water moving into the ocean.
127 
128  ! Local variables
129  real :: Temp_restore ! The temperature that is being restored toward [degC].
130  real :: Salin_restore ! The salinity that is being restored toward [ppt]
131  real :: density_restore ! The potential density that is being restored
132  ! toward [kg m-3].
133  real :: rhoXcp ! The mean density times the heat capacity [J m-3 degC-1].
134  real :: Rho0_mks ! The mean density in MKS units [kg m-3]
135  real :: buoy_rest_const ! A constant relating density anomalies to the
136  ! restoring buoyancy flux [L2 m3 T-3 kg-1 ~> m5 s-3 kg-1].
137 
138  integer :: i, j, is, ie, js, je
139  integer :: isd, ied, jsd, jed
140 
141  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
142  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
143  rho0_mks = cs%Rho0 * us%R_to_kg_m3
144 
145  ! When modifying the code, comment out this error message. It is here
146  ! so that the original (unmodified) version is not accidentally used.
147  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
148  "User forcing routine called without modification." )
149 
150  ! Allocate and zero out the forcing arrays, as necessary. This portion is
151  ! usually not changed.
152  if (cs%use_temperature) then
153  call safe_alloc_ptr(fluxes%evap, isd, ied, jsd, jed)
154  call safe_alloc_ptr(fluxes%lprec, isd, ied, jsd, jed)
155  call safe_alloc_ptr(fluxes%fprec, isd, ied, jsd, jed)
156  call safe_alloc_ptr(fluxes%lrunoff, isd, ied, jsd, jed)
157  call safe_alloc_ptr(fluxes%frunoff, isd, ied, jsd, jed)
158  call safe_alloc_ptr(fluxes%vprec, isd, ied, jsd, jed)
159 
160  call safe_alloc_ptr(fluxes%sw, isd, ied, jsd, jed)
161  call safe_alloc_ptr(fluxes%lw, isd, ied, jsd, jed)
162  call safe_alloc_ptr(fluxes%latent, isd, ied, jsd, jed)
163  call safe_alloc_ptr(fluxes%sens, isd, ied, jsd, jed)
164  else ! This is the buoyancy only mode.
165  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
166  endif
167 
168 
169  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
170 
171  if ( cs%use_temperature ) then
172  ! Set whichever fluxes are to be used here. Any fluxes that
173  ! are always zero do not need to be changed here.
174  do j=js,je ; do i=is,ie
175  ! Fluxes of fresh water through the surface are in units of [R Z T-1 ~> kg m-2 s-1]
176  ! and are positive downward - i.e. evaporation should be negative.
177  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
178  fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
179 
180  ! vprec will be set later, if it is needed for salinity restoring.
181  fluxes%vprec(i,j) = 0.0
182 
183  ! Heat fluxes are in units of W m-2 and are positive into the ocean.
184  fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
185  fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
186  fluxes%sens(i,j) = 0.0 * g%mask2dT(i,j)
187  fluxes%sw(i,j) = 0.0 * g%mask2dT(i,j)
188  enddo ; enddo
189  else ! This is the buoyancy only mode.
190  do j=js,je ; do i=is,ie
191  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
192  ! buoyancy flux is of the same sign as heating the ocean.
193  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
194  enddo ; enddo
195  endif
196 
197  if (cs%restorebuoy) then
198  if (cs%use_temperature) then
199  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
200  ! When modifying the code, comment out this error message. It is here
201  ! so that the original (unmodified) version is not accidentally used.
202  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
203  "Temperature and salinity restoring used without modification." )
204 
205  rhoxcp = rho0_mks * fluxes%C_p
206  do j=js,je ; do i=is,ie
207  ! Set Temp_restore and Salin_restore to the temperature (in degC) and
208  ! salinity (in PSU or ppt) that are being restored toward.
209  temp_restore = 0.0
210  salin_restore = 0.0
211 
212  fluxes%heat_added(i,j) = (g%mask2dT(i,j) * (rhoxcp * us%Z_to_m*us%s_to_T*cs%Flux_const)) * &
213  (temp_restore - sfc_state%SST(i,j))
214  fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
215  ((salin_restore - sfc_state%SSS(i,j)) / (0.5 * (salin_restore + sfc_state%SSS(i,j))))
216  enddo ; enddo
217  else
218  ! When modifying the code, comment out this error message. It is here
219  ! so that the original (unmodified) version is not accidentally used.
220  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
221  "Buoyancy restoring used without modification." )
222 
223  ! The -1 is because density has the opposite sign to buoyancy.
224  buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / rho0_mks
225  do j=js,je ; do i=is,ie
226  ! Set density_restore to an expression for the surface potential
227  ! density [kg m-3] that is being restored toward.
228  density_restore = 1030.0
229 
230  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
231  (density_restore - sfc_state%sfc_density(i,j))
232  enddo ; enddo
233  endif
234  endif ! end RESTOREBUOY
235 

References mom_error_handler::mom_error().

Here is the call graph for this function:

◆ user_surface_forcing_init()

subroutine, public user_surface_forcing::user_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(in), target  diag,
type(user_surface_forcing_cs), pointer  CS 
)

This subroutine initializes the USER_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]diagA structure that is used to regulate diagnostic output.
csA pointer that is set to point to the control structure for this module

Definition at line 240 of file user_surface_forcing.F90.

240  type(time_type), intent(in) :: Time !< The current model time
241  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
242  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
243  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
244  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate diagnostic output.
245  type(user_surface_forcing_CS), pointer :: CS !< A pointer that is set to point to
246  !! the control structure for this module
247 
248 ! This include declares and sets the variable "version".
249 #include "version_variable.h"
250  character(len=40) :: mdl = "user_surface_forcing" ! This module's name.
251 
252  if (associated(cs)) then
253  call mom_error(warning, "USER_surface_forcing_init called with an associated "// &
254  "control structure.")
255  return
256  endif
257  allocate(cs)
258  cs%diag => diag
259 
260  ! Read all relevant parameters and write them to the model log.
261  call log_version(param_file, mdl, version, "")
262  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
263  "If true, Temperature and salinity are used as state "//&
264  "variables.", default=.true.)
265 
266  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
267  "The gravitational acceleration of the Earth.", &
268  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
269  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
270  "The mean ocean density used with BOUSSINESQ true to "//&
271  "calculate accelerations and the mass for conservation "//&
272  "properties, or with BOUSSINSEQ false to convert some "//&
273  "parameters from vertical units of m to kg m-2.", &
274  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
275  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
276  "The background gustiness in the winds.", &
277  units="Pa", default=0.02, scale=us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z)
278 
279  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
280  "If true, the buoyancy fluxes drive the model back "//&
281  "toward some specified surface state with a rate "//&
282  "given by FLUXCONST.", default= .false.)
283  if (cs%restorebuoy) then
284  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
285  "The constant that relates the restoring surface fluxes "//&
286  "to the relative surface anomalies (akin to a piston "//&
287  "velocity). Note the non-MKS units.", &
288  units="m day-1", scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
289  ! Convert CS%Flux_const from m day-1 to m s-1.
290  cs%Flux_const = cs%Flux_const / 86400.0
291  endif
292 

References mom_error_handler::mom_error().

Referenced by mom_surface_forcing::surface_forcing_init().

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

◆ user_wind_forcing()

subroutine, public user_surface_forcing::user_wind_forcing ( 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(user_surface_forcing_cs), pointer  CS 
)

This subroutine sets the surface wind stresses, forcestaux and forcestauy, in [R Z L T-2 ~> Pa]. These are the stresses in the direction of the model grid (i.e. the same direction as the u- and v- velocities).

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
csA pointer to the control structure returned by a previous call to user_surface_forcing_init

Definition at line 52 of file user_surface_forcing.F90.

52  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
53  !! describe the surface state of the ocean.
54  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
55  type(time_type), intent(in) :: day !< The time of the fluxes
56  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
57  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
58  type(user_surface_forcing_CS), pointer :: CS !< A pointer to the control structure returned
59  !! by a previous call to user_surface_forcing_init
60 
61  ! Local variables
62  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
63 
64  ! When modifying the code, comment out this error message. It is here
65  ! so that the original (unmodified) version is not accidentally used.
66  call mom_error(fatal, "User_wind_surface_forcing: " // &
67  "User forcing routine called without modification." )
68 
69  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
70  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
71 
72  ! Allocate the forcing arrays, if necessary.
73  call allocate_mech_forcing(g, forces, stress=.true., ustar=.true.)
74 
75  ! Set the surface wind stresses, in units of [R L Z T-1 ~> Pa]. A positive taux
76  ! accelerates the ocean to the (pseudo-)east.
77 
78  ! The i-loop extends to is-1 so that taux can be used later in the
79  ! calculation of ustar - otherwise the lower bound would be Isq.
80  do j=js,je ; do i=is-1,ieq
81  ! Change this to the desired expression.
82  forces%taux(i,j) = g%mask2dCu(i,j) * 0.0*us%kg_m3_to_R*us%m_s_to_L_T**2*us%L_to_Z
83  enddo ; enddo
84  do j=js-1,jeq ; do i=is,ie
85  forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0 ! Change this to the desired expression.
86  enddo ; enddo
87 
88  ! Set the surface friction velocity, in units of [Z T-1 ~> m s-1]. ustar
89  ! is always positive.
90  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
91  ! This expression can be changed if desired, but need not be.
92  forces%ustar(i,j) = g%mask2dT(i,j) * sqrt((cs%gust_const + &
93  sqrt(0.5*(forces%taux(i-1,j)**2 + forces%taux(i,j)**2) + &
94  0.5*(forces%tauy(i,j-1)**2 + forces%tauy(i,j)**2))) * (us%L_to_Z/cs%Rho0))
95  enddo ; enddo ; endif
96 

References mom_forcing_type::allocate_mech_forcing(), and mom_error_handler::mom_error().

Here is the call graph for this function: