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 [R ~> kg m-3].
133  real :: rhoXcp ! The mean density times the heat capacity [Q R degC-1 ~> J m-3 degC-1].
134  real :: buoy_rest_const ! A constant relating density anomalies to the
135  ! restoring buoyancy flux [L2 m3 T-3 kg-1 ~> m5 s-3 kg-1].
136 
137  integer :: i, j, is, ie, js, je
138  integer :: isd, ied, jsd, jed
139 
140  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
141  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
142 
143  ! When modifying the code, comment out this error message. It is here
144  ! so that the original (unmodified) version is not accidentally used.
145  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
146  "User forcing routine called without modification." )
147 
148  ! Allocate and zero out the forcing arrays, as necessary. This portion is
149  ! usually not changed.
150  if (cs%use_temperature) then
151  call safe_alloc_ptr(fluxes%evap, isd, ied, jsd, jed)
152  call safe_alloc_ptr(fluxes%lprec, isd, ied, jsd, jed)
153  call safe_alloc_ptr(fluxes%fprec, isd, ied, jsd, jed)
154  call safe_alloc_ptr(fluxes%lrunoff, isd, ied, jsd, jed)
155  call safe_alloc_ptr(fluxes%frunoff, isd, ied, jsd, jed)
156  call safe_alloc_ptr(fluxes%vprec, isd, ied, jsd, jed)
157 
158  call safe_alloc_ptr(fluxes%sw, isd, ied, jsd, jed)
159  call safe_alloc_ptr(fluxes%lw, isd, ied, jsd, jed)
160  call safe_alloc_ptr(fluxes%latent, isd, ied, jsd, jed)
161  call safe_alloc_ptr(fluxes%sens, isd, ied, jsd, jed)
162  else ! This is the buoyancy only mode.
163  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
164  endif
165 
166 
167  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
168 
169  if ( cs%use_temperature ) then
170  ! Set whichever fluxes are to be used here. Any fluxes that
171  ! are always zero do not need to be changed here.
172  do j=js,je ; do i=is,ie
173  ! Fluxes of fresh water through the surface are in units of [R Z T-1 ~> kg m-2 s-1]
174  ! and are positive downward - i.e. evaporation should be negative.
175  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
176  fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
177 
178  ! vprec will be set later, if it is needed for salinity restoring.
179  fluxes%vprec(i,j) = 0.0
180 
181  ! Heat fluxes are in units of [Q R Z T-1 ~> W m-2] and are positive into the ocean.
182  fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
183  fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
184  fluxes%sens(i,j) = 0.0 * g%mask2dT(i,j)
185  fluxes%sw(i,j) = 0.0 * g%mask2dT(i,j)
186  enddo ; enddo
187  else ! This is the buoyancy only mode.
188  do j=js,je ; do i=is,ie
189  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
190  ! buoyancy flux is of the same sign as heating the ocean.
191  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
192  enddo ; enddo
193  endif
194 
195  if (cs%restorebuoy) then
196  if (cs%use_temperature) then
197  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
198  ! When modifying the code, comment out this error message. It is here
199  ! so that the original (unmodified) version is not accidentally used.
200  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
201  "Temperature and salinity restoring used without modification." )
202 
203  rhoxcp = cs%Rho0 * fluxes%C_p
204  do j=js,je ; do i=is,ie
205  ! Set Temp_restore and Salin_restore to the temperature (in degC) and
206  ! salinity (in PSU or ppt) that are being restored toward.
207  temp_restore = 0.0
208  salin_restore = 0.0
209 
210  fluxes%heat_added(i,j) = (g%mask2dT(i,j) * (rhoxcp * cs%Flux_const)) * &
211  (temp_restore - sfc_state%SST(i,j))
212  fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
213  ((salin_restore - sfc_state%SSS(i,j)) / (0.5 * (salin_restore + sfc_state%SSS(i,j))))
214  enddo ; enddo
215  else
216  ! When modifying the code, comment out this error message. It is here
217  ! so that the original (unmodified) version is not accidentally used.
218  call mom_error(fatal, "User_buoyancy_surface_forcing: " // &
219  "Buoyancy restoring used without modification." )
220 
221  ! The -1 is because density has the opposite sign to buoyancy.
222  buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / cs%Rho0
223  do j=js,je ; do i=is,ie
224  ! Set density_restore to an expression for the surface potential
225  ! density [kg m-3] that is being restored toward.
226  density_restore = 1030.0*us%kg_m3_to_R
227 
228  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
229  (density_restore - sfc_state%sfc_density(i,j))
230  enddo ; enddo
231  endif
232  endif ! end RESTOREBUOY
233 

◆ 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 238 of file user_surface_forcing.F90.

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

◆ 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