MOM6
dumbbell_surface_forcing.F90
Go to the documentation of this file.
1 !> Surface forcing for the dumbbell test case
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
8 use mom_domains, only : pass_var, pass_vector, agrid
9 use mom_error_handler, only : mom_error, fatal, warning, is_root_pe
12 use mom_grid, only : ocean_grid_type
13 use mom_io, only : file_exists, read_data
15 use mom_time_manager, only : time_type, operator(+), operator(/), get_time
19 use mom_variables, only : surface
20 
21 implicit none ; private
22 
24 
25 !> Control structure for the dumbbell test case forcing
26 type, public :: dumbbell_surface_forcing_cs ; private
27  logical :: use_temperature !< If true, temperature and salinity are used as
28  !! state variables.
29  logical :: restorebuoy !< If true, use restoring surface buoyancy forcing.
30  real :: rho0 !< The density used in the Boussinesq approximation [R ~> kg m-3].
31  real :: g_earth !< The gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
32  real :: flux_const !< The restoring rate at the surface [Z T-1 ~> m s-1].
33  real :: gust_const !< A constant unresolved background gustiness
34  !! that contributes to ustar [Pa].
35  real :: slp_amplitude !< The amplitude of pressure loading [Pa] applied
36  !! to the reservoirs
37  real :: slp_period !< Period of sinusoidal pressure wave
38  real, dimension(:,:), allocatable :: &
39  forcing_mask !< A mask regulating where forcing occurs
40  real, dimension(:,:), allocatable :: &
41  s_restore !< The surface salinity field toward which to
42  !! restore [ppt].
43  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to regulate the
44  !! timing of diagnostic output.
46 
47 contains
48 
49 !> Surface buoyancy (heat and fresh water) fluxes for the dumbbell test case
50 subroutine dumbbell_buoyancy_forcing(state, fluxes, day, dt, G, US, CS)
51  type(surface), intent(inout) :: state !< A structure containing fields that
52  !! describe the surface state of the ocean.
53  type(forcing), intent(inout) :: fluxes !< A structure containing pointers to any
54  !! possible forcing fields. Unused fields
55  !! have NULL ptrs.
56  type(time_type), intent(in) :: day !< Time of the fluxes.
57  real, intent(in) :: dt !< The amount of time over which
58  !! the fluxes apply [s]
59  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
60  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
61  type(dumbbell_surface_forcing_cs), pointer :: cs !< A control structure returned by a previous
62  !! call to dumbbell_surface_forcing_init
63  ! Local variables
64  real :: temp_restore ! The temperature that is being restored toward [degC].
65  real :: salin_restore ! The salinity that is being restored toward [ppt].
66  real :: density_restore ! The potential density that is being restored
67  ! toward [kg m-3].
68  real :: rhoxcp ! The mean density times the heat capacity [J m-3 degC-1].
69  integer :: i, j, is, ie, js, je
70  integer :: isd, ied, jsd, jed
71 
72  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
73  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
74 
75 
76  ! Allocate and zero out the forcing arrays, as necessary.
77  if (cs%use_temperature) then
78  call safe_alloc_ptr(fluxes%evap, isd, ied, jsd, jed)
79  call safe_alloc_ptr(fluxes%lprec, isd, ied, jsd, jed)
80  call safe_alloc_ptr(fluxes%fprec, isd, ied, jsd, jed)
81  call safe_alloc_ptr(fluxes%lrunoff, isd, ied, jsd, jed)
82  call safe_alloc_ptr(fluxes%frunoff, isd, ied, jsd, jed)
83  call safe_alloc_ptr(fluxes%vprec, isd, ied, jsd, jed)
84 
85  call safe_alloc_ptr(fluxes%sw, isd, ied, jsd, jed)
86  call safe_alloc_ptr(fluxes%lw, isd, ied, jsd, jed)
87  call safe_alloc_ptr(fluxes%latent, isd, ied, jsd, jed)
88  call safe_alloc_ptr(fluxes%sens, isd, ied, jsd, jed)
89  else ! This is the buoyancy only mode.
90  call safe_alloc_ptr(fluxes%buoy, isd, ied, jsd, jed)
91  endif
92 
93 
94  ! MODIFY THE CODE IN THE FOLLOWING LOOPS TO SET THE BUOYANCY FORCING TERMS.
95 
96  if ( cs%use_temperature ) then
97  ! Set whichever fluxes are to be used here. Any fluxes that
98  ! are always zero do not need to be changed here.
99  do j=js,je ; do i=is,ie
100  ! Fluxes of fresh water through the surface are in units of [kg m-2 s-1]
101  ! and are positive downward - i.e. evaporation should be negative.
102  fluxes%evap(i,j) = -0.0 * g%mask2dT(i,j)
103  fluxes%lprec(i,j) = 0.0 * g%mask2dT(i,j)
104 
105  ! vprec will be set later, if it is needed for salinity restoring.
106  fluxes%vprec(i,j) = 0.0
107 
108  ! Heat fluxes are in units of [W m-2] and are positive into the ocean.
109  fluxes%lw(i,j) = 0.0 * g%mask2dT(i,j)
110  fluxes%latent(i,j) = 0.0 * g%mask2dT(i,j)
111  fluxes%sens(i,j) = 0.0 * g%mask2dT(i,j)
112  fluxes%sw(i,j) = 0.0 * g%mask2dT(i,j)
113  enddo ; enddo
114  else ! This is the buoyancy only mode.
115  do j=js,je ; do i=is,ie
116  ! fluxes%buoy is the buoyancy flux into the ocean [L2 T-3 ~> m2 s-3]. A positive
117  ! buoyancy flux is of the same sign as heating the ocean.
118  fluxes%buoy(i,j) = 0.0 * g%mask2dT(i,j)
119  enddo ; enddo
120  endif
121 
122  if (cs%use_temperature .and. cs%restorebuoy) then
123  do j=js,je ; do i=is,ie
124  ! Set density_restore to an expression for the surface potential
125  ! density [kg m-3] that is being restored toward.
126  if (cs%forcing_mask(i,j)>0.) then
127  fluxes%vprec(i,j) = - (g%mask2dT(i,j) * (cs%Rho0*cs%Flux_const)) * &
128  ((cs%S_restore(i,j) - state%SSS(i,j)) / (0.5 * (cs%S_restore(i,j) + state%SSS(i,j))))
129 
130  endif
131  enddo ; enddo
132  endif
133 
134 end subroutine dumbbell_buoyancy_forcing
135 
136 !> Dynamic forcing for the dumbbell test case
137 subroutine dumbbell_dynamic_forcing(state, fluxes, day, dt, G, CS)
138  type(surface), intent(inout) :: state !< A structure containing fields that
139  !! describe the surface state of the ocean.
140  type(forcing), intent(inout) :: fluxes !< A structure containing pointers to any
141  !! possible forcing fields. Unused fields
142  !! have NULL ptrs.
143  type(time_type), intent(in) :: day !< Time of the fluxes.
144  real, intent(in) :: dt !< The amount of time over which
145  !! the fluxes apply [s]
146  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
147  type(dumbbell_surface_forcing_cs), pointer :: cs !< A control structure returned by a previous
148  !! call to dumbbell_surface_forcing_init
149  ! Local variables
150  integer :: i, j, is, ie, js, je
151  integer :: isd, ied, jsd, jed
152  integer :: idays, isecs
153  real :: deg_rad, rdays
154 
155 
156  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
157  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
158 
159  deg_rad = atan(1.0)*4.0/180.
160 
161  call get_time(day,isecs,idays)
162  rdays = real(idays) + real(isecs)/8.64e4
163  ! This could be: rdays = time_type_to_real(day)/8.64e4
164 
165  ! Allocate and zero out the forcing arrays, as necessary.
166  call safe_alloc_ptr(fluxes%p_surf, isd, ied, jsd, jed)
167  call safe_alloc_ptr(fluxes%p_surf_full, isd, ied, jsd, jed)
168 
169  do j=js,je ; do i=is,ie
170  fluxes%p_surf(i,j) = cs%forcing_mask(i,j)* cs%slp_amplitude * &
171  g%mask2dT(i,j) * sin(deg_rad*(rdays/cs%slp_period))
172  fluxes%p_surf_full(i,j) = cs%forcing_mask(i,j) * cs%slp_amplitude * &
173  g%mask2dT(i,j) * sin(deg_rad*(rdays/cs%slp_period))
174  enddo ; enddo
175 
176 end subroutine dumbbell_dynamic_forcing
177 
178 !> Reads and sets up the forcing for the dumbbell test case
179 subroutine dumbbell_surface_forcing_init(Time, G, US, param_file, diag, CS)
180  type(time_type), intent(in) :: time !< The current model time.
181  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
182  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
183  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
184  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to
185  !! regulate diagnostic output.
187  pointer :: cs !< A pointer to the control structure for this module
188  ! Local variables
189  real :: s_surf, s_range
190  real :: x, y
191  integer :: i, j
192 #include "version_variable.h"
193  character(len=40) :: mdl = "dumbbell_surface_forcing" ! This module's name.
194 
195  if (associated(cs)) then
196  call mom_error(warning, "dumbbell_surface_forcing_init called with an associated "// &
197  "control structure.")
198  return
199  endif
200  allocate(cs)
201  cs%diag => diag
202 
203  ! Read all relevant parameters and write them to the model log.
204  call log_version(param_file, mdl, version, "")
205  call get_param(param_file, mdl, "ENABLE_THERMODYNAMICS", cs%use_temperature, &
206  "If true, Temperature and salinity are used as state "//&
207  "variables.", default=.true.)
208 
209  call get_param(param_file, mdl, "G_EARTH", cs%G_Earth, &
210  "The gravitational acceleration of the Earth.", &
211  units="m s-2", default = 9.80, scale=us%m_to_L**2*us%Z_to_m*us%T_to_s**2)
212  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
213  "The mean ocean density used with BOUSSINESQ true to "//&
214  "calculate accelerations and the mass for conservation "//&
215  "properties, or with BOUSSINSEQ false to convert some "//&
216  "parameters from vertical units of m to kg m-2.", &
217  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
218  call get_param(param_file, mdl, "DUMBBELL_SLP_AMP", cs%slp_amplitude, &
219  "Amplitude of SLP forcing in reservoirs.", &
220  units="kg m2 s-1", default = 10000.0)
221  call get_param(param_file, mdl, "DUMBBELL_SLP_PERIOD", cs%slp_period, &
222  "Periodicity of SLP forcing in reservoirs.", &
223  units="days", default = 1.0)
224  call get_param(param_file, mdl, "DUMBBELL_SLP_PERIOD", cs%slp_period, &
225  "Periodicity of SLP forcing in reservoirs.", &
226  units="days", default = 1.0)
227  call get_param(param_file, mdl,"INITIAL_SSS", s_surf, &
228  "Initial surface salinity", units="1e-3", default=34.0, do_not_log=.true.)
229  call get_param(param_file, mdl,"INITIAL_S_RANGE", s_range, &
230  "Initial salinity range (bottom - surface)", units="1e-3", &
231  default=2., do_not_log=.true.)
232 
233  call get_param(param_file, mdl, "RESTOREBUOY", cs%restorebuoy, &
234  "If true, the buoyancy fluxes drive the model back "//&
235  "toward some specified surface state with a rate "//&
236  "given by FLUXCONST.", default= .false.)
237  if (cs%restorebuoy) then
238  call get_param(param_file, mdl, "FLUXCONST", cs%Flux_const, &
239  "The constant that relates the restoring surface fluxes "//&
240  "to the relative surface anomalies (akin to a piston "//&
241  "velocity). Note the non-MKS units.", &
242  units="m day-1", scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
243  ! Convert CS%Flux_const from m day-1 to m s-1.
244  cs%Flux_const = cs%Flux_const / 86400.0
245 
246 
247  allocate(cs%forcing_mask(g%isd:g%ied, g%jsd:g%jed)); cs%forcing_mask(:,:)=0.0
248  allocate(cs%S_restore(g%isd:g%ied, g%jsd:g%jed))
249 
250  do j=g%jsc,g%jec
251  do i=g%isc,g%iec
252  ! Compute normalized zonal coordinates (x,y=0 at center of domain)
253  x = ( g%geoLonT(i,j) - g%west_lon ) / g%len_lon - 0.5
254  y = ( g%geoLatT(i,j) - g%south_lat ) / g%len_lat - 0.5
255  cs%forcing_mask(i,j)=0
256  cs%S_restore(i,j) = s_surf
257  if ((x>0.25)) then
258  cs%forcing_mask(i,j) = 1
259  cs%S_restore(i,j) = s_surf + s_range
260  elseif ((x<-0.25)) then
261  cs%forcing_mask(i,j) = 1
262  cs%S_restore(i,j) = s_surf - s_range
263  endif
264  enddo
265  enddo
266  endif
267 
268 end subroutine dumbbell_surface_forcing_init
269 
270 end module dumbbell_surface_forcing
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_forcing_type::allocate_forcing_type
subroutine, public allocate_forcing_type(G, fluxes, water, heat, ustar, press, shelf, iceberg, salt)
Conditionally allocate fields within the forcing type.
Definition: MOM_forcing_type.F90:2811
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_safe_alloc
Convenience functions for safely allocating memory without accidentally reallocating pointer and caus...
Definition: MOM_safe_alloc.F90:3
mom_diag_mediator::query_averaging_enabled
logical function, public query_averaging_enabled(diag_cs, time_int, time_end)
Call this subroutine to determine whether the averaging is currently enabled. .true....
Definition: MOM_diag_mediator.F90:1850
dumbbell_surface_forcing::dumbbell_surface_forcing_init
subroutine, public dumbbell_surface_forcing_init(Time, G, US, param_file, diag, CS)
Reads and sets up the forcing for the dumbbell test case.
Definition: dumbbell_surface_forcing.F90:180
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_diag_mediator::register_diag_field
integer function, public register_diag_field(module_name, field_name, axes_in, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, cell_methods, x_cell_method, y_cell_method, v_cell_method, conversion, v_extensive)
Returns the "diag_mediator" handle for a group (native, CMOR, z-coord, ...) of diagnostics derived fr...
Definition: MOM_diag_mediator.F90:1878
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
dumbbell_surface_forcing::dumbbell_buoyancy_forcing
subroutine, public dumbbell_buoyancy_forcing(state, fluxes, day, dt, G, US, CS)
Surface buoyancy (heat and fresh water) fluxes for the dumbbell test case.
Definition: dumbbell_surface_forcing.F90:51
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_tracer_flow_control
Orchestrates the registration and calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:2
dumbbell_surface_forcing
Surface forcing for the dumbbell test case.
Definition: dumbbell_surface_forcing.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_tracer_flow_control::tracer_flow_control_cs
The control structure for orchestrating the calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:75
mom_error_handler::is_root_pe
logical function, public is_root_pe()
This returns .true. if the current PE is the root PE.
Definition: MOM_error_handler.F90:44
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:50
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
dumbbell_surface_forcing::dumbbell_surface_forcing_cs
Control structure for the dumbbell test case forcing.
Definition: dumbbell_surface_forcing.F90:26
mom_error_handler::mom_error
subroutine, public mom_error(level, message, all_print)
This provides a convenient interface for writing an mpp_error message with run-time filter based on a...
Definition: MOM_error_handler.F90:72
mom_safe_alloc::safe_alloc_ptr
Allocate a pointer to a 1-d, 2-d or 3-d array.
Definition: MOM_safe_alloc.F90:12
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
dumbbell_surface_forcing::dumbbell_dynamic_forcing
subroutine, public dumbbell_dynamic_forcing(state, fluxes, day, dt, G, CS)
Dynamic forcing for the dumbbell test case.
Definition: dumbbell_surface_forcing.F90:138
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
mom_tracer_flow_control::call_tracer_set_forcing
subroutine, public call_tracer_set_forcing(state, fluxes, day_start, day_interval, G, CS)
This subroutine calls the individual tracer modules' subroutines to specify or read quantities relate...
Definition: MOM_tracer_flow_control.F90:382
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239