MOM6
neverland_surface_forcing Module Reference

Detailed Description

Wind and buoyancy forcing for the Neverland configurations.

Data Types

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

Functions/Subroutines

subroutine, public neverland_wind_forcing (sfc_state, forces, day, G, US, CS)
 Sets the surface wind stresses, forcestaux and forcestauy for the Neverland forcing configuration. More...
 
real function cosbell (x, L)
 Returns the value of a cosine-bell function evaluated at x/L. More...
 
real function spike (x, L)
 Returns the value of a sin-spike function evaluated at x/L. More...
 
subroutine, public neverland_buoyancy_forcing (sfc_state, fluxes, day, dt, G, US, CS)
 Surface fluxes of buoyancy for the Neverland configurations. More...
 
subroutine, public neverland_surface_forcing_init (Time, G, US, param_file, diag, CS)
 Initializes the Neverland control structure. More...
 

Function/Subroutine Documentation

◆ cosbell()

real function neverland_surface_forcing::cosbell ( real, intent(in)  x,
real, intent(in)  L 
)
private

Returns the value of a cosine-bell function evaluated at x/L.

Parameters
[in]xnon-dimensional position
[in]lnon-dimensional width

Definition at line 116 of file Neverland_surface_forcing.F90.

116 
117  real , intent(in) :: x !< non-dimensional position
118  real , intent(in) :: L !< non-dimensional width
119  real :: PI !< 3.1415926... calculated as 4*atan(1)
120 
121  pi = 4.0*atan(1.0)
122  cosbell = 0.5 * (1 + cos(pi*min(abs(x/l),1.0)))

◆ neverland_buoyancy_forcing()

subroutine, public neverland_surface_forcing::neverland_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(inout)  G,
type(unit_scale_type), intent(in)  US,
type(neverland_surface_forcing_cs), pointer  CS 
)

Surface fluxes of buoyancy for the Neverland configurations.

Parameters
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in,out]fluxesForcing fields.
[in]dayTime used for determining the fluxes.
[in]dtForcing time step (s).
[in,out]gGrid structure.
[in]usA dimensional unit scaling type
csControl structure for this module.

Definition at line 139 of file Neverland_surface_forcing.F90.

139  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
140  !! describe the surface state of the ocean.
141  type(forcing), intent(inout) :: fluxes !< Forcing fields.
142  type(time_type), intent(in) :: day !< Time used for determining the fluxes.
143  real, intent(in) :: dt !< Forcing time step (s).
144  type(ocean_grid_type), intent(inout) :: G !< Grid structure.
145  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
146  type(Neverland_surface_forcing_CS), pointer :: CS !< Control structure for this module.
147  ! Local variables
148  real :: buoy_rest_const ! A constant relating density anomalies to the
149  ! restoring buoyancy flux [L2 m3 T-3 kg-1 ~> m5 s-3 kg-1].
150  real :: density_restore ! De
151  integer :: i, j, is, ie, js, je
152  integer :: isd, ied, jsd, jed
153 
154  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
155  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
156 
157  ! When modifying the code, comment out this error message. It is here
158  ! so that the original (unmodified) version is not accidentally used.
159 
160  ! Allocate and zero out the forcing arrays, as necessary. This portion is
161  ! usually not changed.
162  if (cs%use_temperature) then
163  call mom_error(fatal, "Neverland_buoyancy_forcing: " // &
164  "Temperature and salinity mode not coded!" )
165  else
166  ! This is the buoyancy only mode.
167  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
168  endif
169 
170 
171  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
172  if (cs%restorebuoy .and. cs%first_call) then
173  call safe_alloc_ptr(cs%buoy_restore, isd, ied, jsd, jed)
174  cs%first_call = .false.
175  ! Set CS%buoy_restore(i,j) here
176  endif
177 
178  if ( cs%use_temperature ) then
179  call mom_error(fatal, "Neverland_buoyancy_surface_forcing: " // &
180  "Temperature/salinity restoring not coded!" )
181  else ! This is the buoyancy only mode.
182  do j=js,je ; do i=is,ie
183  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
184  ! buoyancy flux is of the same sign as heating the ocean.
185  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
186  enddo ; enddo
187  endif
188 
189  if (cs%restorebuoy) then
190  if (cs%use_temperature) then
191  call mom_error(fatal, "Neverland_buoyancy_surface_forcing: " // &
192  "Temperature/salinity restoring not coded!" )
193  else
194  ! When modifying the code, comment out this error message. It is here
195  ! so that the original (unmodified) version is not accidentally used.
196 
197  ! The -1 is because density has the opposite sign to buoyancy.
198  buoy_rest_const = -1.0 * (cs%G_Earth * us%m_to_Z*us%T_to_s*cs%Flux_const) / cs%Rho0
199  do j=js,je ; do i=is,ie
200  ! Set density_restore to an expression for the surface potential
201  ! density [kg m-3] that is being restored toward.
202  density_restore = 1030.0
203 
204  fluxes%buoy(i,j) = g%mask2dT(i,j) * buoy_rest_const * &
205  (density_restore - sfc_state%sfc_density(i,j))
206  enddo ; enddo
207  endif
208  endif ! end RESTOREBUOY
209 

◆ neverland_surface_forcing_init()

subroutine, public neverland_surface_forcing::neverland_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(neverland_surface_forcing_cs), pointer  CS 
)

Initializes the Neverland control structure.

Parameters
[in]timeThe current model time.
[in]gThe ocean's grid structure.
[in]usA dimensional unit scaling type
[in]param_fileA structure indicating the open file to parse for model parameter values.
[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 214 of file Neverland_surface_forcing.F90.

214  type(time_type), intent(in) :: Time !< The current model time.
215  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
216  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
217  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to parse for
218  !! model parameter values.
219  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate diagnostic output.
220  type(Neverland_surface_forcing_CS), pointer :: CS !< A pointer that is set to point to the control structure
221  !! for this module
222  ! This include declares and sets the variable "version".
223 #include "version_variable.h"
224  ! Local variables
225  character(len=40) :: mdl = "Neverland_surface_forcing" ! This module's name.
226 
227  if (associated(cs)) then
228  call mom_error(warning, "Neverland_surface_forcing_init called with an associated "// &
229  "control structure.")
230  return
231  endif
232  allocate(cs)
233  cs%diag => diag
234 
235  ! Read all relevant parameters and write them to the model log.
236  call log_version(param_file, mdl, version, "")
237  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
238  "If true, Temperature and salinity are used as state "//&
239  "variables.", default=.true.)
240 
241  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
242  "The gravitational acceleration of the Earth.", &
243  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
244  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
245  "The mean ocean density used with BOUSSINESQ true to "//&
246  "calculate accelerations and the mass for conservation "//&
247  "properties, or with BOUSSINSEQ false to convert some "//&
248  "parameters from vertical units of m to kg m-2.", &
249  units="kg m-3", default=1035.0)
250 ! call get_param(param_file, mdl, "GUST_CONST", CS%gust_const, &
251 ! "The background gustiness in the winds.", units="Pa", &
252 ! default=0.02)
253 
254  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
255  "If true, the buoyancy fluxes drive the model back "//&
256  "toward some specified surface state with a rate "//&
257  "given by FLUXCONST.", default= .false.)
258 
259  if (cs%restorebuoy) then
260  call get_param(param_file, mdl, "FLUXCONST", cs%flux_const, &
261  "The constant that relates the restoring surface fluxes "//&
262  "to the relative surface anomalies (akin to a piston "//&
263  "velocity). Note the non-MKS units.", units="m day-1", &
264  fail_if_missing=.true.)
265  ! Convert CS%flux_const from m day-1 to m s-1.
266  cs%flux_const = cs%flux_const / 86400.0
267  endif
268 

◆ neverland_wind_forcing()

subroutine, public neverland_surface_forcing::neverland_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(neverland_surface_forcing_cs), pointer  CS 
)

Sets the surface wind stresses, forcestaux and forcestauy for the Neverland forcing configuration.

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]dayTime used for determining the fluxes.
[in,out]gGrid structure.
[in]usA dimensional unit scaling type
csControl structure for this module.

Definition at line 51 of file Neverland_surface_forcing.F90.

51  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
52  !! describe the surface state of the ocean.
53  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
54  type(time_type), intent(in) :: day !< Time used for determining the fluxes.
55  type(ocean_grid_type), intent(inout) :: G !< Grid structure.
56  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
57  type(Neverland_surface_forcing_CS), pointer :: CS !< Control structure for this module.
58 
59  ! Local variables
60  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
61  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
62  real :: x, y
63  real :: PI
64  real :: tau_max, off
65 
66  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
67  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
68  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
69  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
70 
71  ! Allocate the forcing arrays, if necessary.
72  call allocate_mech_forcing(g, forces, stress=.true.)
73 
74  ! Set the surface wind stresses, in units of Pa. A positive taux
75  ! accelerates the ocean to the (pseudo-)east.
76 
77  ! The i-loop extends to is-1 so that taux can be used later in the
78  ! calculation of ustar - otherwise the lower bound would be Isq.
79  pi = 4.0*atan(1.0)
80  forces%taux(:,:) = 0.0
81  tau_max = 0.2
82  off = 0.02
83  do j=js,je ; do i=is-1,ieq
84 ! x = (G%geoLonT(i,j)-G%west_lon)/G%len_lon
85  y = (g%geoLatT(i,j)-g%south_lat)/g%len_lat
86 ! forces%taux(I,j) = G%mask2dCu(I,j) * 0.0
87 
88  if (y <= 0.29) then
89  forces%taux(i,j) = forces%taux(i,j) + tau_max * ( (1/0.29)*y - ( 1/(2*pi) )*sin( (2*pi*y) / 0.29 ) )
90  endif
91  if ((y > 0.29) .and. (y <= (0.8-off))) then
92  forces%taux(i,j) = forces%taux(i,j) + tau_max *(0.35+0.65*cos(pi*(y-0.29)/(0.51-off)) )
93  endif
94  if ((y > (0.8-off)) .and. (y <= (1-off))) then
95  forces%taux(i,j) = forces%taux(i,j) + tau_max *( 1.5*( (y-1+off) - (0.1/pi)*sin(10.0*pi*(y-0.8+off)) ) )
96  endif
97  enddo ; enddo
98 
99  do j=js-1,jeq ; do i=is,ie
100  forces%tauy(i,j) = g%mask2dCv(i,j) * 0.0 ! Change this to the desired expression.
101  enddo ; enddo
102 
103  ! Set the surface friction velocity, in units of m s-1. ustar
104  ! is always positive.
105 ! if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
106 ! ! This expression can be changed if desired, but need not be.
107 ! forces%ustar(i,j) = US%m_to_Z*US%T_to_s * G%mask2dT(i,j) * sqrt(CS%gust_const/CS%Rho0 + &
108 ! sqrt(0.5*(forces%taux(I-1,j)**2 + forces%taux(I,j)**2) + &
109 ! 0.5*(forces%tauy(i,J-1)**2 + forces%tauy(i,J)**2))/CS%Rho0)
110 ! enddo ; enddo ; endif
111 

◆ spike()

real function neverland_surface_forcing::spike ( real, intent(in)  x,
real, intent(in)  L 
)
private

Returns the value of a sin-spike function evaluated at x/L.

Parameters
[in]xnon-dimensional position
[in]lnon-dimensional width

Definition at line 127 of file Neverland_surface_forcing.F90.

127 
128  real , intent(in) :: x !< non-dimensional position
129  real , intent(in) :: L !< non-dimensional width
130  real :: PI !< 3.1415926... calculated as 4*atan(1)
131 
132  pi = 4.0*atan(1.0)
133  spike = (1 - sin(pi*min(abs(x/l),0.5)))