MOM6
meso_surface_forcing Module Reference

Detailed Description

Sets forcing for the MESO configuration.

Rewritten by Robert Hallberg, June 2009.

This file contains the subroutines that a user should modify to to set the surface wind stresses and fluxes of buoyancy or temperature and fresh water. They are called when the run-time parameters WIND_CONFIG or BUOY_CONFIG are set to "USER". The standard version has simple examples, along with run-time error messages that will cause the model to abort if this code has not been modified. This code is intended for use with relatively simple specifications of the forcing. For more complicated forms, it is probably a good idea to read the forcing from input files using "file" for WIND_CONFIG and BUOY_CONFIG.

MESO_buoyancy forcing is used to set the surface buoyancy forcing, which may include a number of fresh water flux fields (evap, liq_precip, froz_precip, liq_runoff, froz_runoff, and vprec) and the surface heat fluxes (sw, lw, latent and sens) if temperature and salinity are state variables, or it may simply be the buoyancy flux if it is not. This routine also has coded a restoring to surface values of temperature and salinity.

Data Types

type  meso_surface_forcing_cs
 This control structure is used to store parameters associated with the MESO forcing. More...
 

Functions/Subroutines

subroutine, public meso_buoyancy_forcing (sfc_state, fluxes, day, dt, G, US, CS)
 This subroutine sets up the MESO buoyancy forcing, which uses control-theory style specification restorative buoyancy fluxes at large scales. More...
 
subroutine, public meso_surface_forcing_init (Time, G, US, param_file, diag, CS)
 Initialize the MESO surface forcing module. More...
 

Variables

logical first_call = .true.
 True until after the first call to the MESO forcing routines.
 

Function/Subroutine Documentation

◆ meso_buoyancy_forcing()

subroutine, public meso_surface_forcing::meso_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(meso_surface_forcing_cs), pointer  CS 
)

This subroutine sets up the MESO buoyancy forcing, which uses control-theory style specification restorative buoyancy fluxes at large scales.

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 MESO_surface_forcing_init

Definition at line 59 of file MESO_surface_forcing.F90.

59  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
60  !! describe the surface state of the ocean.
61  type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
62  type(time_type), intent(in) :: day !< The time of the fluxes
63  real, intent(in) :: dt !< The amount of time over which
64  !! the fluxes apply [s]
65  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
66  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
67  type(MESO_surface_forcing_CS), pointer :: CS !< A pointer to the control structure returned by
68  !! a previous call to MESO_surface_forcing_init
69 
70 ! When temperature is used, there are long list of fluxes that need to be
71 ! set - essentially the same as for a full coupled model, but most of these
72 ! can be simply set to zero. The net fresh water flux should probably be
73 ! set in fluxes%evap and fluxes%lprec, with any salinity restoring
74 ! appearing in fluxes%vprec, and the other water flux components
75 ! (fprec, lrunoff and frunoff) left as arrays full of zeros.
76 ! Evap is usually negative and precip is usually positive. All heat fluxes
77 ! are in W m-2 and positive for heat going into the ocean. All fresh water
78 ! fluxes are in kg m-2 s-1 and positive for water moving into the ocean.
79 
80  real :: Temp_restore ! The temperature that is being restored toward [degC].
81  real :: Salin_restore ! The salinity that is being restored toward [ppt]
82  real :: density_restore ! The potential density that is being restored toward [R ~> kg m-3].
83  real :: rhoXcp ! The mean density times the heat capacity [Q R degC-1 ~> J m-3 degC-1].
84  real :: buoy_rest_const ! A constant relating density anomalies to the
85  ! restoring buoyancy flux [L2 T-3 R-1 ~> m5 s-3 kg-1].
86 
87  integer :: i, j, is, ie, js, je
88  integer :: isd, ied, jsd, jed
89 
90  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
91  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
92 
93  ! When modifying the code, comment out this error message. It is here
94  ! so that the original (unmodified) version is not accidentally used.
95 
96  ! Allocate and zero out the forcing arrays, as necessary. This portion is
97  ! usually not changed.
98  if (cs%use_temperature) then
99  call safe_alloc_ptr(fluxes%evap, isd, ied, jsd, jed)
100  call safe_alloc_ptr(fluxes%lprec, isd, ied, jsd, jed)
101  call safe_alloc_ptr(fluxes%fprec, isd, ied, jsd, jed)
102  call safe_alloc_ptr(fluxes%lrunoff, isd, ied, jsd, jed)
103  call safe_alloc_ptr(fluxes%frunoff, isd, ied, jsd, jed)
104  call safe_alloc_ptr(fluxes%vprec, isd, ied, jsd, jed)
105 
106  call safe_alloc_ptr(fluxes%sw, isd, ied, jsd, jed)
107  call safe_alloc_ptr(fluxes%lw, isd, ied, jsd, jed)
108  call safe_alloc_ptr(fluxes%latent, isd, ied, jsd, jed)
109  call safe_alloc_ptr(fluxes%sens, isd, ied, jsd, jed)
110  call safe_alloc_ptr(fluxes%heat_content_lprec, isd, ied, jsd, jed)
111  else ! This is the buoyancy only mode.
112  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
113  endif
114 
115 
116  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
117  if (cs%restorebuoy .and. first_call) then !#CTRL# .or. associated(CS%ctrl_forcing_CSp)) then
118  call safe_alloc_ptr(cs%T_Restore, isd, ied, jsd, jed)
119  call safe_alloc_ptr(cs%S_Restore, isd, ied, jsd, jed)
120  call safe_alloc_ptr(cs%Heat, isd, ied, jsd, jed)
121  call safe_alloc_ptr(cs%PmE, isd, ied, jsd, jed)
122  call safe_alloc_ptr(cs%Solar, isd, ied, jsd, jed)
123 
124  call mom_read_data(trim(cs%inputdir)//trim(cs%SSTrestore_file), "SST", &
125  cs%T_Restore(:,:), g%Domain)
126  call mom_read_data(trim(cs%inputdir)//trim(cs%salinityrestore_file), "SAL", &
127  cs%S_Restore(:,:), g%Domain)
128  call mom_read_data(trim(cs%inputdir)//trim(cs%heating_file), "Heat", &
129  cs%Heat(:,:), g%Domain, scale=us%W_m2_to_QRZ_T)
130  call mom_read_data(trim(cs%inputdir)//trim(cs%PmE_file), "PmE", &
131  cs%PmE(:,:), g%Domain, scale=us%m_to_Z*us%T_to_s)
132  call mom_read_data(trim(cs%inputdir)//trim(cs%Solar_file), "NET_SOL", &
133  cs%Solar(:,:), g%Domain, scale=us%W_m2_to_QRZ_T)
134  first_call = .false.
135  endif
136 
137  if ( cs%use_temperature ) then
138  ! Set whichever fluxes are to be used here. Any fluxes that
139  ! are always zero do not need to be changed here.
140  do j=js,je ; do i=is,ie
141  ! Fluxes of fresh water through the surface are in units of [kg m-2 s-1]
142  ! and are positive downward - i.e. evaporation should be negative.
143  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
144  fluxes%lprec(i,j) = cs%PmE(i,j) * cs%Rho0 * g%mask2dT(i,j)
145 
146  ! vprec will be set later, if it is needed for salinity restoring.
147  fluxes%vprec(i,j) = 0.0
148 
149  ! Heat fluxes are in units of [Q R Z T-1 ~> W m-2] and are positive into the ocean.
150  fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
151  fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
152  fluxes%sens(i,j) = cs%Heat(i,j) * g%mask2dT(i,j)
153  fluxes%sw(i,j) = cs%Solar(i,j) * g%mask2dT(i,j)
154  enddo ; enddo
155  else ! This is the buoyancy only mode.
156  do j=js,je ; do i=is,ie
157  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
158  ! buoyancy flux is of the same sign as heating the ocean.
159  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
160  enddo ; enddo
161  endif
162 
163  if (cs%restorebuoy) then
164  if (cs%use_temperature) then
165  call safe_alloc_ptr(fluxes%heat_added, isd, ied, jsd, jed)
166  ! When modifying the code, comment out this error message. It is here
167  ! so that the original (unmodified) version is not accidentally used.
168 ! call MOM_error(FATAL, "MESO_buoyancy_surface_forcing: " // &
169 ! "Temperature and salinity restoring used without modification." )
170 
171  rhoxcp = cs%Rho0 * fluxes%C_p
172  do j=js,je ; do i=is,ie
173  ! Set Temp_restore and Salin_restore to the temperature (in degC) and
174  ! salinity (in ppt or PSU) that are being restored toward.
175  if (g%mask2dT(i,j) > 0) then
176  fluxes%heat_added(i,j) = g%mask2dT(i,j) * &
177  ((cs%T_Restore(i,j) - sfc_state%SST(i,j)) * rhoxcp * cs%Flux_const)
178  fluxes%vprec(i,j) = - (cs%Rho0 * cs%Flux_const) * &
179  (cs%S_Restore(i,j) - sfc_state%SSS(i,j)) / &
180  (0.5*(sfc_state%SSS(i,j) + cs%S_Restore(i,j)))
181  else
182  fluxes%heat_added(i,j) = 0.0
183  fluxes%vprec(i,j) = 0.0
184  endif
185  enddo ; enddo
186  else
187  ! When modifying the code, comment out this error message. It is here
188  ! so that the original (unmodified) version is not accidentally used.
189  call mom_error(fatal, "MESO_buoyancy_surface_forcing: " // &
190  "Buoyancy restoring used without modification." )
191 
192  ! The -1 is because density has the opposite sign to buoyancy.
193  buoy_rest_const = -1.0 * (cs%G_Earth * cs%Flux_const) / cs%Rho0
194  do j=js,je ; do i=is,ie
195  ! Set density_restore to an expression for the surface potential
196  ! density [R ~> kg m-3] that is being restored toward.
197  density_restore = 1030.0 * us%kg_m3_to_R
198 
199  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
200  (density_restore - sfc_state%sfc_density(i,j))
201  enddo ; enddo
202  endif
203  endif ! end RESTOREBUOY
204 

◆ meso_surface_forcing_init()

subroutine, public meso_surface_forcing::meso_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(meso_surface_forcing_cs), pointer  CS 
)

Initialize the MESO 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
csA pointer that is set to point to the control structure for this module

Definition at line 209 of file MESO_surface_forcing.F90.

209 
210  type(time_type), intent(in) :: Time !< The current model time
211  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
212  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
213  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
214  type(diag_ctrl), target, intent(inout) :: diag !< structure used to regulate diagnostic output
215  type(MESO_surface_forcing_CS), pointer :: CS !< A pointer that is set to point to the
216  !! control structure for this module
217 
218 ! This include declares and sets the variable "version".
219 #include "version_variable.h"
220  character(len=40) :: mdl = "MESO_surface_forcing" ! This module's name.
221 
222  if (associated(cs)) then
223  call mom_error(warning, "MESO_surface_forcing_init called with an associated "// &
224  "control structure.")
225  return
226  endif
227  allocate(cs)
228  cs%diag => diag
229 
230  ! Read all relevant parameters and write them to the model log.
231  call log_version(param_file, mdl, version, "")
232  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
233  "If true, Temperature and salinity are used as state "//&
234  "variables.", default=.true.)
235 
236  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
237  "The gravitational acceleration of the Earth.", &
238  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
239  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
240  "The mean ocean density used with BOUSSINESQ true to "//&
241  "calculate accelerations and the mass for conservation "//&
242  "properties, or with BOUSSINSEQ false to convert some "//&
243  "parameters from vertical units of m to kg m-2.", &
244  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
245  call get_param(param_file, mdl, "GUST_CONST", cs%gust_const, &
246  "The background gustiness in the winds.", units="Pa", default=0.0)
247 
248  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
249  "If true, the buoyancy fluxes drive the model back "//&
250  "toward some specified surface state with a rate "//&
251  "given by FLUXCONST.", default= .false.)
252 
253  if (cs%restorebuoy) then
254  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
255  "The constant that relates the restoring surface fluxes to the relative "//&
256  "surface anomalies (akin to a piston velocity). Note the non-MKS units.", &
257  default=0.0, units="m day-1", scale=us%m_to_Z/(86400.0*us%s_to_T))
258 
259  call get_param(param_file, mdl, "SSTRESTORE_FILE", cs%SSTrestore_file, &
260  "The file with the SST toward which to restore in "//&
261  "variable TEMP.", fail_if_missing=.true.)
262  call get_param(param_file, mdl, "SALINITYRESTORE_FILE", cs%salinityrestore_file, &
263  "The file with the surface salinity toward which to "//&
264  "restore in variable SALT.", fail_if_missing=.true.)
265  call get_param(param_file, mdl, "SENSIBLEHEAT_FILE", cs%heating_file, &
266  "The file with the non-shortwave heat flux in "//&
267  "variable Heat.", fail_if_missing=.true.)
268  call get_param(param_file, mdl, "PRECIP_FILE", cs%PmE_file, &
269  "The file with the net precipiation minus evaporation "//&
270  "in variable PmE.", fail_if_missing=.true.)
271  call get_param(param_file, mdl, "SHORTWAVE_FILE", cs%Solar_file, &
272  "The file with the shortwave heat flux in "//&
273  "variable NET_SOL.", fail_if_missing=.true.)
274  call get_param(param_file, mdl, "INPUTDIR", cs%inputdir, default=".")
275  cs%inputdir = slasher(cs%inputdir)
276 
277  endif
278