MOM_barotropic.F90

1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> Barotropic solver
7
8use mom_checksums, only : chksum0
9use mom_coms, only : any_across_pes
10use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end, clock_routine
11use mom_debugging, only : hchksum, uvchksum, bchksum
14use mom_domains, only : min_across_pes, clone_mom_domain, deallocate_mom_domain
15use mom_domains, only : to_all, scalar_pair, agrid, corner, mom_domain_type
16use mom_domains, only : create_group_pass, do_group_pass, group_pass_type
17use mom_domains, only : start_group_pass, complete_group_pass, pass_var, pass_vector
18use mom_error_handler, only : mom_error, mom_mesg, fatal, warning, is_root_pe
19use mom_file_parser, only : get_param, log_param, log_version, param_file_type
21use mom_grid, only : ocean_grid_type
24use mom_io, only : vardesc, var_desc, mom_read_data, slasher, north_face, east_face
26use mom_open_boundary, only : obc_direction_e, obc_direction_w
27use mom_open_boundary, only : obc_direction_n, obc_direction_s, obc_segment_type
28use mom_restart, only : register_restart_field, register_restart_pair
29use mom_restart, only : query_initialized, mom_restart_cs
31use mom_self_attr_load, only : sal_cs
33use mom_time_manager, only : time_type, real_to_time, operator(+), operator(-)
39
40implicit none ; private
41
42#include <MOM_memory.h>
43#ifdef STATIC_MEMORY_
44# ifndef BTHALO_
45# define BTHALO_ 0
46# endif
47# define WHALOI_ MAX(BTHALO_-NIHALO_,0)
48# define WHALOJ_ MAX(BTHALO_-NJHALO_,0)
49# define NIMEMW_ 1-WHALOI_:NIMEM_+WHALOI_
50# define NJMEMW_ 1-WHALOJ_:NJMEM_+WHALOJ_
51# define NIMEMBW_ -WHALOI_:NIMEM_+WHALOI_
52# define NJMEMBW_ -WHALOJ_:NJMEM_+WHALOJ_
53# define SZIW_(G) NIMEMW_
54# define SZJW_(G) NJMEMW_
55# define SZIBW_(G) NIMEMBW_
56# define SZJBW_(G) NJMEMBW_
57#else
58# define NIMEMW_ :
59# define NJMEMW_ :
60# define NIMEMBW_ :
61# define NJMEMBW_ :
62# define SZIW_(G) G%isdw:G%iedw
63# define SZJW_(G) G%jsdw:G%jedw
64# define SZIBW_(G) G%isdw-1:G%iedw
65# define SZJBW_(G) G%jsdw-1:G%jedw
66#endif
67
70
71! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
72! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
73! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
74! vary with the Boussinesq approximation, the Boussinesq variant is given first.
75
76!> The barotropic stepping open boundary condition type
77type, private :: bt_obc_type
78 real, allocatable :: cg_u(:,:) !< The external wave speed at u-points [L T-1 ~> m s-1].
79 real, allocatable :: cg_v(:,:) !< The external wave speed at u-points [L T-1 ~> m s-1].
80 real, allocatable :: dz_u(:,:) !< The total vertical column extent at the u-points [Z ~> m].
81 real, allocatable :: dz_v(:,:) !< The total vertical column extent at the v-points [Z ~> m].
82 real, allocatable :: uhbt(:,:) !< The zonal barotropic thickness fluxes specified
83 !! for open boundary conditions (if any) [H L2 T-1 ~> m3 s-1 or kg s-1].
84 real, allocatable :: vhbt(:,:) !< The meridional barotropic thickness fluxes specified
85 !! for open boundary conditions (if any) [H L2 T-1 ~> m3 s-1 or kg s-1].
86 real, allocatable :: ubt_outer(:,:) !< The zonal velocities just outside the domain,
87 !! as set by the open boundary conditions [L T-1 ~> m s-1].
88 real, allocatable :: vbt_outer(:,:) !< The meridional velocities just outside the domain,
89 !! as set by the open boundary conditions [L T-1 ~> m s-1].
90 real, allocatable :: ssh_outer_u(:,:) !< The surface height outside of the domain
91 !! at a u-point with an open boundary condition [Z ~> m].
92 real, allocatable :: ssh_outer_v(:,:) !< The surface height outside of the domain
93 !! at a v-point with an open boundary condition [Z ~> m].
94 integer, allocatable :: u_obc_type(:,:) !< An integer encoding the type and direction of u-point OBCs
95 integer, allocatable :: v_obc_type(:,:) !< An integer encoding the type and direction of v-point OBCs
96 logical :: u_obcs_on_pe !< True if this PE has an open boundary at any u-points.
97 logical :: v_obcs_on_pe !< True if this PE has an open boundary at any v-points.
98 !>@{ Index ranges on the local PE for the open boundary conditions in various directions
99 integer :: is_u_w_obc, ie_u_w_obc, js_u_w_obc, je_u_w_obc
100 integer :: is_u_e_obc, ie_u_e_obc, js_u_e_obc, je_u_e_obc
101 integer :: is_v_s_obc, ie_v_s_obc, js_v_s_obc, je_v_s_obc
102 integer :: is_v_n_obc, ie_v_n_obc, js_v_n_obc, je_v_n_obc
103 !>@}
104
105 type(group_pass_type) :: pass_uv !< Structure for group halo pass of vectors
106 type(group_pass_type) :: scalar_pass !< Structure for group halo pass of scalars
107end type bt_obc_type
108
109integer, parameter :: specified_obc = 1 !< An integer used to encode a specified OBC point
110integer, parameter :: flather_obc = 2 !< An integer used to encode a Flather OBC point
111integer, parameter :: gradient_obc = 4 !< An integer used to encode a gradient OBC point
112
113!> The barotropic stepping control structure
114type, public :: barotropic_cs ; private
115 real allocable_, dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: frhatu
116 !< The fraction of the total column thickness interpolated to u grid points in each layer [nondim].
117 real allocable_, dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: frhatv
118 !< The fraction of the total column thickness interpolated to v grid points in each layer [nondim].
119 real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: idatu
120 !< Inverse of the total thickness at u grid points [H-1 ~> m-1 or m2 kg-1].
121 real, allocatable, dimension(:,:) :: lin_drag_u
122 !< A spatially varying linear drag coefficient acting on the zonal barotropic flow
123 !! [H T-1 ~> m s-1 or kg m-2 s-1].
124 real, allocatable, dimension(:,:) :: ubt_ic
125 !< The barotropic solvers estimate of the zonal velocity that will be the initial
126 !! condition for the next call to btstep [L T-1 ~> m s-1].
127 real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: ubtav
128 !< The barotropic zonal velocity averaged over the baroclinic time step [L T-1 ~> m s-1].
129 real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: idatv
130 !< Inverse of the basin depth at v grid points [Z-1 ~> m-1].
131 real, allocatable, dimension(:,:) :: lin_drag_v
132 !< A spatially varying linear drag coefficient acting on the zonal barotropic flow
133 !! [H T-1 ~> m s-1 or kg m-2 s-1].
134 real, allocatable, dimension(:,:) :: vbt_ic
135 !< The barotropic solvers estimate of the zonal velocity that will be the initial
136 !! condition for the next call to btstep [L T-1 ~> m s-1].
137 real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: vbtav
138 !< The barotropic meridional velocity averaged over the baroclinic time step [L T-1 ~> m s-1].
139 real allocable_, dimension(NIMEM_,NJMEM_) :: eta_cor
140 !< The difference between the free surface height from the barotropic calculation and the sum
141 !! of the layer thicknesses. This difference is imposed as a forcing term in the barotropic
142 !! calculation over a baroclinic timestep [H ~> m or kg m-2].
143 real, allocatable, dimension(:,:) :: eta_cor_bound
144 !< A limit on the rate at which eta_cor can be applied while avoiding instability
145 !! [H T-1 ~> m s-1 or kg m-2 s-1]. This is only used if CS%bound_BT_corr is true.
146 real allocable_, dimension(NIMEMW_,NJMEMW_) :: &
147 ua_polarity, & !< Test vector components for checking grid polarity [nondim]
148 va_polarity, & !< Test vector components for checking grid polarity [nondim]
149 bathyt !< A copy of bathyT (ocean bottom depth) with wide halos [Z ~> m]
150 real allocable_, dimension(NIMEMW_,NJMEMW_) :: iareat
151 !< This is a copy of G%IareaT with wide halos, but will
152 !! still utilize the macro IareaT when referenced, [L-2 ~> m-2].
153 real allocable_, dimension(NIMEMBW_,NJMEMW_) :: &
154 dy_cu, & !< A copy of G%dy_Cu with wide halos [L ~> m].
155 idxcu, & !< A copy of G%IdxCu with wide halos [L-1 ~> m-1].
156 obcmask_u !< An array to multiplicatively mask out changes at OBC points, 0 or 1 [nondim]
157 real allocable_, dimension(NIMEMW_,NJMEMBW_) :: &
158 dx_cv, & !< A copy of G%dx_Cv with wide halos [L ~> m].
159 idycv, & !< A copy of G%IdyCv with wide halos [L-1 ~> m-1].
160 obcmask_v !< An array to multiplicatively mask out changes at OBC points, 0 or 1 [nondim]
161 real, allocatable, dimension(:,:) :: &
162 d_u_cor, & !< A simply averaged depth at u points recast as a thickness [H ~> m or kg m-2]
163 d_v_cor, & !< A simply averaged depth at v points recast as a thickness [H ~> m or kg m-2]
164 q_d !< f / D at PV points [T-1 H-1 ~> s-1 m-1 or m2 s-1 kg-1]
165 real, allocatable, dimension(:,:,:) :: &
166 q_wt !< The area weights for the thicknesses around a corner point to be used when
167 !! calculating PV for use in the Coriolis term, taking OBCs into account [L2 ~> m2].
168 !! The order of the 4 values at a point is the order in which the neighboring
169 !! tracer points occur in memory, i.e. SW, SE, NW then NE.
170 real, allocatable :: frhatu1(:,:,:) !< Predictor step values of frhatu stored for diagnostics [nondim]
171 real, allocatable :: frhatv1(:,:,:) !< Predictor step values of frhatv stored for diagnostics [nondim]
172 real, allocatable :: iareat_obcmask(:,:) !< If non-zero, work on given points [L-2 ~> m-2].
173
174 type(bt_obc_type) :: bt_obc !< A structure with all of this modules fields
175 !! for applying open boundary conditions.
176
177 real :: dtbt !< The barotropic time step [T ~> s].
178 real :: dtbt_fraction !< The fraction of the maximum time-step that
179 !! should used [nondim]. The default is 0.98.
180 real :: dtbt_max !< The maximum stable barotropic time step [T ~> s].
181 real :: dt_bt_filter !< The time-scale over which the barotropic mode solutions are
182 !! filtered [T ~> s] if positive, or as a fraction of DT if
183 !! negative [nondim]. This can never be taken to be longer than 2*dt.
184 !! Set this to 0 to apply no filtering.
185 integer :: nstep_last = 0 !< The number of barotropic timesteps per baroclinic
186 !! time step the last time btstep was called.
187 real :: bebt !< A nondimensional number, from 0 to 1, that
188 !! determines the gravity wave time stepping scheme [nondim].
189 !! 0.0 gives a forward-backward scheme, while 1.0
190 !! give backward Euler. In practice, bebt should be
191 !! of order 0.2 or greater.
192 real :: rho_bt_lin !< A density that is used to convert total water column thicknesses
193 !! into mass in non-Boussinesq mode with linearized options in the
194 !! barotropic solver or when estimating the stable barotropic timestep
195 !! without access to the full baroclinic model state [R ~> kg m-3]
196 logical :: split !< If true, use the split time stepping scheme.
197 logical :: bound_bt_corr !< If true, the magnitude of the fake mass source
198 !! in the barotropic equation that drives the two
199 !! estimates of the free surface height toward each
200 !! other is bounded to avoid driving corrective
201 !! velocities that exceed MAXCFL_BT_CONT.
202 logical :: gradual_bt_ics !< If true, adjust the initial conditions for the
203 !! barotropic solver to the values from the layered
204 !! solution over a whole timestep instead of
205 !! instantly. This is a decent approximation to the
206 !! inclusion of sum(u dh_dt) while also correcting
207 !! for truncation errors.
208 logical :: sadourny !< If true, the Coriolis terms are discretized
209 !! with Sadourny's energy conserving scheme,
210 !! otherwise the Arakawa & Hsu scheme is used. If
211 !! the deformation radius is not resolved Sadourny's
212 !! scheme should probably be used.
213 logical :: integral_bt_cont !< If true, use the time-integrated velocity over the barotropic steps
214 !! to determine the integrated transports used to update the continuity
215 !! equation. Otherwise the transports are the sum of the transports
216 !! based on a series of instantaneous velocities and the BT_CONT_TYPE
217 !! for transports. This is only valid if a BT_CONT_TYPE is used.
218 logical :: bt_adjust_src_for_filter !< If true, increases the rate at which BT mass sources are
219 !! applied so that they are all used up before the steps within the
220 !! filtering period start. This avoids the mass sink driving the SSH
221 !! below the bottom during the period of filtering.
222 logical :: bt_limit_integral_transport !< If true, limit the time-integrated transports by the
223 !! initial volume accounting for sinks of mass.
224 logical :: integral_obcs !< This is true if integral_bt_cont is true and there are open boundary
225 !! conditions being applied somewhere in the global domain.
226 logical :: nonlinear_continuity !< If true, the barotropic continuity equation
227 !! uses the full ocean thickness for transport.
228 integer :: nonlin_cont_update_period !< The number of barotropic time steps
229 !! between updates to the face area, or 0 only to
230 !! update at the start of a call to btstep. The
231 !! default is 1.
232 logical :: bt_project_velocity !< If true, step the barotropic velocity first
233 !! and project out the velocity tendency by 1+BEBT
234 !! when calculating the transport. The default
235 !! (false) is to use a predictor continuity step to
236 !! find the pressure field, and then do a corrector
237 !! continuity step using a weighted average of the
238 !! old and new velocities, with weights of (1-BEBT) and BEBT.
239 logical :: nonlin_stress !< If true, use the full depth of the ocean at the start of the
240 !! barotropic step when calculating the surface stress contribution to
241 !! the barotropic accelerations. Otherwise use the depth based on bathyT.
242 real :: bt_coriolis_scale !< A factor by which the barotropic Coriolis acceleration anomaly
243 !! terms are scaled [nondim].
244 integer :: answer_date !< The vintage of the expressions in the barotropic solver.
245 !! Values below 20190101 recover the answers from the end of 2018,
246 !! while higher values use more efficient or general expressions.
247
248 logical :: dynamic_psurf !< If true, add a dynamic pressure due to a viscous
249 !! ice shelf, for instance.
250 real :: dmin_dyn_psurf !< The minimum total thickness to use in limiting the size
251 !! of the dynamic surface pressure for stability [H ~> m or kg m-2].
252 real :: ice_strength_length !< The length scale at which the damping rate
253 !! due to the ice strength should be the same as if
254 !! a Laplacian were applied [L ~> m].
255 real :: const_dyn_psurf !< The constant that scales the dynamic surface
256 !! pressure [nondim]. Stable values are < ~1.0.
257 !! The default is 0.9.
258 logical :: calculate_sal !< If true, calculate self-attraction and loading.
259 logical :: tidal_sal_bug !< If true, the tidal self-attraction and loading anomaly in the
260 !! barotropic solver has the wrong sign, replicating a long-standing
261 !! bug.
262 real :: g_extra !< A nondimensional factor by which gtot is enhanced [nondim].
263 integer :: hvel_scheme !< An integer indicating how the thicknesses at
264 !! velocity points are calculated. Valid values are
265 !! given by the parameters defined below:
266 !! HARMONIC, ARITHMETIC, HYBRID, and FROM_BT_CONT
267 logical :: strong_drag !< If true, use a stronger estimate of the retarding
268 !! effects of strong bottom drag.
269 logical :: rescale_strong_drag !< If true, reduce the barotropic contribution to the layer
270 !! accelerations to account for the difference between the forces that
271 !! can be counteracted by the stronger drag with BT_STRONG_DRAG and the
272 !! average of the layer viscous remnants after a baroclinic timestep.
273 logical :: linear_wave_drag !< If true, apply a linear drag to the barotropic
274 !! velocities, using rates set by lin_drag_u & _v
275 !! divided by the depth of the ocean.
276 logical :: linearized_bt_pv !< If true, the PV and interface thicknesses used
277 !! in the barotropic Coriolis calculation is time
278 !! invariant and linearized.
279 logical :: use_filter !< If true, use streaming band-pass filter to detect the
280 !! instantaneous tidal signals in the simulation.
281 logical :: linear_freq_drag !< If true, apply a linear frequency-dependent drag to the tidal
282 !! velocities. The streaming band-pass filter must be turned on.
283 logical :: use_wide_halos !< If true, use wide halos and march in during the
284 !! barotropic time stepping for efficiency.
285 integer :: min_stencil !< The minimum stencil width to use with the wide halo iterations.
286 !! A nonzero value may reflect the distribution of OBC faces or it
287 !! may be useful for debugging purposes.
288 logical :: clip_velocity !< If true, limit any velocity components that are
289 !! are large enough for a CFL number to exceed
290 !! CFL_trunc. This should only be used as a
291 !! desperate debugging measure.
292 logical :: debug !< If true, write verbose checksums for debugging purposes.
293 logical :: debug_bt !< If true, write verbose checksums from within the barotropic
294 !! time-stepping loop for debugging purposes.
295 logical :: debug_wide_halos !< If true, write the checksums on the full wide halos. Otherwise
296 !! only the output for the final computational domain is written.
297 real :: vel_underflow !< Velocity components smaller than vel_underflow
298 !! are set to 0 [L T-1 ~> m s-1].
299 real :: maxvel !< Velocity components greater than maxvel are
300 !! truncated to maxvel [L T-1 ~> m s-1].
301 real :: cfl_trunc !< If clip_velocity is true, velocity components will
302 !! be truncated when they are large enough that the
303 !! corresponding CFL number exceeds this value [nondim].
304 real :: maxcfl_bt_cont !< The maximum permitted CFL number associated with the
305 !! barotropic accelerations from the summed velocities
306 !! times the time-derivatives of thicknesses [nondim]. The
307 !! default is 0.1, and there will probably be real
308 !! problems if this were set close to 1.
309 logical :: bt_cont_bounds !< If true, use the BT_cont_type variables to set limits
310 !! on the magnitude of the corrective mass fluxes.
311 logical :: visc_rem_u_uh0 !< If true, use the viscous remnants when estimating
312 !! the barotropic velocities that were used to
313 !! calculate uh0 and vh0. False is probably the
314 !! better choice.
315 logical :: adjust_bt_cont !< If true, adjust the curve fit to the BT_cont type
316 !! that is used by the barotropic solver to match the
317 !! transport about which the flow is being linearized.
318 logical :: use_old_coriolis_bracket_bug !< If True, use an order of operations
319 !! that is not bitwise rotationally symmetric in the
320 !! meridional Coriolis term of the barotropic solver.
321 logical :: tidal_sal_flather !< Apply adjustment to external gravity wave speed
322 !! consistent with tidal self-attraction and loading
323 !! used within the barotropic solver
324 logical :: wt_uv_bug = .true. !< If true, recover a bug that wt_[uv] that is not normalized.
325 logical :: exterior_obc_bug = .true. !< If true, recover a bug with boundary conditions
326 !! inside the domain.
327 logical :: interior_obc_pv !< If true, use only interior ocean points at OBCs to specify the PV
328 !! used in the barotropic Coriolis anomalies. Otherwise the
329 !! calculation relies on bathymetry and eta being projected outward
330 !! across OBCs. Unfortunately, this option does change answers near
331 !! convex (peninsula-type) pairs of OBC segments.
332 type(time_type), pointer :: time => null() !< A pointer to the ocean models clock.
333 type(diag_ctrl), pointer :: diag => null() !< A structure that is used to regulate
334 !! the timing of diagnostic output.
335 type(mom_domain_type), pointer :: bt_domain => null() !< Barotropic MOM domain
336 type(hor_index_type), pointer :: debug_bt_hi => null() !< debugging copy of horizontal index_type
337 type(sal_cs), pointer :: sal_csp => null() !< Control structure for SAL
338 type(harmonic_analysis_cs), pointer :: ha_csp => null() !< Control structure for harmonic analysis
339 type(filter_cs) :: filt_cs_u, & !< Control structures for the streaming band-pass filter of ubt
340 filt_cs_v !< Control structures for the streaming band-pass filter of vbt
341 type(wave_drag_cs) :: drag_cs !< Control structures for the frequency-dependent drag
342 logical :: module_is_initialized = .false. !< If true, module has been initialized
343
344 integer :: isdw !< The lower i-memory limit for the wide halo arrays.
345 integer :: iedw !< The upper i-memory limit for the wide halo arrays.
346 integer :: jsdw !< The lower j-memory limit for the wide halo arrays.
347 integer :: jedw !< The upper j-memory limit for the wide halo arrays.
348
349 type(group_pass_type) :: pass_q_dcor !< Handle for a group halo pass
350 type(group_pass_type) :: pass_gtot !< Handle for a group halo pass
351 type(group_pass_type) :: pass_tmp_uv !< Handle for a group halo pass
352 type(group_pass_type) :: pass_eta_bt_rem !< Handle for a group halo pass
353 type(group_pass_type) :: pass_force_hbt0_cor_ref !< Handle for a group halo pass
354 type(group_pass_type) :: pass_dat_uv !< Handle for a group halo pass
355 type(group_pass_type) :: pass_eta_ubt !< Handle for a group halo pass
356 type(group_pass_type) :: pass_etaav !< Handle for a group halo pass
357 type(group_pass_type) :: pass_ubt_cor !< Handle for a group halo pass
358 type(group_pass_type) :: pass_ubta_uhbta !< Handle for a group halo pass
359 type(group_pass_type) :: pass_e_anom !< Handle for a group halo pass
360 type(group_pass_type) :: pass_spv_avg !< Handle for a group halo pass
361
362 !>@{ Diagnostic IDs
363 integer :: id_pfu_bt = -1, id_pfv_bt = -1, id_coru_bt = -1, id_corv_bt = -1
364 integer :: id_ldu_bt = -1, id_ldv_bt = -1, id_eta_cor = -1
365 integer :: id_ubtforce = -1, id_vbtforce = -1, id_uaccel = -1, id_vaccel = -1
366 integer :: id_visc_rem_u = -1, id_visc_rem_v = -1, id_bt_rem_u = -1, id_bt_rem_v = -1
367 integer :: id_ubt = -1, id_vbt = -1, id_eta_bt = -1, id_ubtav = -1, id_vbtav = -1
368 integer :: id_ubt_st = -1, id_vbt_st = -1, id_eta_st = -1
369 integer :: id_ubtdt = -1, id_vbtdt = -1
370 integer :: id_ubt_hifreq = -1, id_vbt_hifreq = -1, id_eta_hifreq = -1
371 integer :: id_uhbt_hifreq = -1, id_vhbt_hifreq = -1, id_eta_pred_hifreq = -1
372 integer :: id_etapf_hifreq = -1, id_etapf_anom = -1
373 integer :: id_gtotn = -1, id_gtots = -1, id_gtote = -1, id_gtotw = -1
374 integer :: id_uhbt = -1, id_frhatu = -1, id_vhbt = -1, id_frhatv = -1
375 integer :: id_frhatu1 = -1, id_frhatv1 = -1
376
377 integer :: id_btc_fa_u_ee = -1, id_btc_fa_u_e0 = -1, id_btc_fa_u_w0 = -1, id_btc_fa_u_ww = -1
378 integer :: id_btc_ubt_ee = -1, id_btc_ubt_ww = -1
379 integer :: id_btc_fa_v_nn = -1, id_btc_fa_v_n0 = -1, id_btc_fa_v_s0 = -1, id_btc_fa_v_ss = -1
380 integer :: id_btc_vbt_nn = -1, id_btc_vbt_ss = -1
381 integer :: id_btc_fa_u_rat0 = -1, id_btc_fa_v_rat0 = -1, id_btc_fa_h_rat0 = -1
382 integer :: id_uhbt0 = -1, id_vhbt0 = -1
383 integer :: id_ssh_u_obc = -1, id_ssh_v_obc = -1, id_ubt_obc = -1, id_vbt_obc = -1
384 !>@}
385
386end type barotropic_cs
387
388!> A description of the functional dependence of transport at a u-point
389type, private :: local_bt_cont_u_type
390 real :: fa_u_ee !< The effective open face area for zonal barotropic transport
391 !! drawing from locations far to the east [H L ~> m2 or kg m-1].
392 real :: fa_u_e0 !< The effective open face area for zonal barotropic transport
393 !! drawing from nearby to the east [H L ~> m2 or kg m-1].
394 real :: fa_u_w0 !< The effective open face area for zonal barotropic transport
395 !! drawing from nearby to the west [H L ~> m2 or kg m-1].
396 real :: fa_u_ww !< The effective open face area for zonal barotropic transport
397 !! drawing from locations far to the west [H L ~> m2 or kg m-1].
398 real :: ubt_ww !< uBT_WW is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY
399 !! the time-integrated barotropic velocity [L ~> m], beyond which the marginal
400 !! open face area is FA_u_WW. uBT_WW must be non-negative.
401 real :: ubt_ee !< uBT_EE is a barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY
402 !! the time-integrated barotropic velocity [L ~> m], beyond which the marginal
403 !! open face area is FA_u_EE. uBT_EE must be non-positive.
404 real :: uh_crvw !< The curvature of face area with velocity for flow from the west [H T2 L-1 ~> s2 or kg s2 m-3]
405 !! or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.
406 real :: uh_crve !< The curvature of face area with velocity for flow from the east [H T2 L-1 ~> s2 or kg s2 m-3]
407 !! or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.
408 real :: uh_ww !< The zonal transport when ubt=ubt_WW [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent
409 !! time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].
410 real :: uh_ee !< The zonal transport when ubt=ubt_EE [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent
411 !! time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].
413
414!> A description of the functional dependence of transport at a v-point
415type, private :: local_bt_cont_v_type
416 real :: fa_v_nn !< The effective open face area for meridional barotropic transport
417 !! drawing from locations far to the north [H L ~> m2 or kg m-1].
418 real :: fa_v_n0 !< The effective open face area for meridional barotropic transport
419 !! drawing from nearby to the north [H L ~> m2 or kg m-1].
420 real :: fa_v_s0 !< The effective open face area for meridional barotropic transport
421 !! drawing from nearby to the south [H L ~> m2 or kg m-1].
422 real :: fa_v_ss !< The effective open face area for meridional barotropic transport
423 !! drawing from locations far to the south [H L ~> m2 or kg m-1].
424 real :: vbt_ss !< vBT_SS is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY
425 !! the time-integrated barotropic velocity [L ~> m], beyond which the marginal
426 !! open face area is FA_v_SS. vBT_SS must be non-negative.
427 real :: vbt_nn !< vBT_NN is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY
428 !! the time-integrated barotropic velocity [L ~> m], beyond which the marginal
429 !! open face area is FA_v_NN. vBT_NN must be non-positive.
430 real :: vh_crvs !< The curvature of face area with velocity for flow from the south [H T2 L-1 ~> s2 or kg s2 m-3]
431 !! or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.
432 real :: vh_crvn !< The curvature of face area with velocity for flow from the north [H T2 L-1 ~> s2 or kg s2 m-3]
433 !! or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.
434 real :: vh_ss !< The meridional transport when vbt=vbt_SS [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent
435 !! time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].
436 real :: vh_nn !< The meridional transport when vbt=vbt_NN [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent
437 !! time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].
439
440!> A container for passing around active tracer point memory limits
441type, private :: memory_size_type
442 !>@{ Currently active memory limits
443 integer :: isdw, iedw, jsdw, jedw ! The memory limits of the wide halo arrays.
444 !>@}
445end type memory_size_type
446
447!>@{ CPU time clock IDs
448integer :: id_clock_sync=-1, id_clock_calc=-1
449integer :: id_clock_calc_pre=-1, id_clock_calc_post=-1
450integer :: id_clock_pass_step=-1, id_clock_pass_pre=-1, id_clock_pass_post=-1
451!>@}
452
453!>@{ Enumeration values for various schemes
454integer, parameter :: harmonic = 1
455integer, parameter :: arithmetic = 2
456integer, parameter :: hybrid = 3
457integer, parameter :: from_bt_cont = 4
458integer, parameter :: hybrid_bt_cont = 5
459character*(20), parameter :: hybrid_string = "HYBRID"
460character*(20), parameter :: harmonic_string = "HARMONIC"
461character*(20), parameter :: arithmetic_string = "ARITHMETIC"
462character*(20), parameter :: bt_cont_string = "FROM_BT_CONT"
463!>@}
464
465!> A negligible parameter which avoids division by zero, but is too small to
466!! modify physical values [nondim].
467real, parameter :: subroundoff = 1e-30
468
469contains
470
471!> This subroutine time steps the barotropic equations explicitly.
472!! For gravity waves, anything between a forwards-backwards scheme
473!! and a simulated backwards Euler scheme is used, with bebt between
474!! 0.0 and 1.0 determining the scheme. In practice, bebt must be of
475!! order 0.2 or greater. A forwards-backwards treatment of the
476!! Coriolis terms is always used.
477subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce, &
478 eta_PF_in, U_Cor, V_Cor, accel_layer_u, accel_layer_v, &
479 eta_out, uhbtav, vhbtav, G, GV, US, CS, &
480 visc_rem_u, visc_rem_v, SpV_avg, ADp, OBC, BT_cont, eta_PF_start, &
481 taux_bot, tauy_bot, uh0, vh0, u_uh0, v_vh0, etaav)
482 type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
483 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
484 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
485 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(in) :: u_in !< The initial (3-D) zonal
486 !! velocity [L T-1 ~> m s-1].
487 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(in) :: v_in !< The initial (3-D) meridional
488 !! velocity [L T-1 ~> m s-1].
489 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: eta_in !< The initial barotropic free surface height
490 !! anomaly or column mass anomaly [H ~> m or kg m-2].
491 real, intent(in) :: dt !< The time increment to integrate over [T ~> s].
492 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(in) :: bc_accel_u !< The zonal baroclinic accelerations,
493 !! [L T-2 ~> m s-2].
494 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(in) :: bc_accel_v !< The meridional baroclinic accelerations,
495 !! [L T-2 ~> m s-2].
496 type(mech_forcing), intent(in) :: forces !< A structure with the driving mechanical forces
497 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: pbce !< The baroclinic pressure anomaly in each layer
498 !! due to free surface height anomalies
499 !! [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
500 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: eta_pf_in !< The 2-D eta field (either SSH anomaly or
501 !! column mass anomaly) that was used to calculate the input
502 !! pressure gradient accelerations (or its final value if
503 !! eta_PF_start is provided [H ~> m or kg m-2].
504 !! Note: eta_in, pbce, and eta_PF_in must have up-to-date
505 !! values in the first point of their halos.
506 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(in) :: u_cor !< The (3-D) zonal velocities used to
507 !! calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1].
508 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(in) :: v_cor !< The (3-D) meridional velocities used to
509 !! calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1].
510 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(out) :: accel_layer_u !< The zonal acceleration of each layer due
511 !! to the barotropic calculation [L T-2 ~> m s-2].
512 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(out) :: accel_layer_v !< The meridional acceleration of each layer
513 !! due to the barotropic calculation [L T-2 ~> m s-2].
514 real, dimension(SZI_(G),SZJ_(G)), intent(out) :: eta_out !< The final barotropic free surface
515 !! height anomaly or column mass anomaly [H ~> m or kg m-2].
516 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: uhbtav !< the barotropic zonal volume or mass
517 !! fluxes averaged through the barotropic steps
518 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
519 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: vhbtav !< the barotropic meridional volume or mass
520 !! fluxes averaged through the barotropic steps
521 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
522 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
523 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(in) :: visc_rem_u !< Both the fraction of the momentum
524 !! originally in a layer that remains after a time-step of
525 !! viscosity, and the fraction of a time-step's worth of a
526 !! barotropic acceleration that a layer experiences after
527 !! viscosity is applied, in the zonal direction [nondim].
528 !! Visc_rem_u is between 0 (at the bottom) and 1 (far above).
529 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(in) :: visc_rem_v !< Ditto for meridional direction [nondim].
530 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: spv_avg !< The column average specific volume, used
531 !! in non-Boussinesq OBC calculations [R-1 ~> m3 kg-1]
532 type(accel_diag_ptrs), pointer :: adp !< Acceleration diagnostic pointers
533 type(ocean_obc_type), pointer :: obc !< The open boundary condition structure.
534 type(bt_cont_type), pointer :: bt_cont !< A structure with elements that describe
535 !! the effective open face areas as a function of barotropic
536 !! flow.
537 real, dimension(:,:), pointer :: eta_pf_start !< The eta field consistent with the pressure
538 !! gradient at the start of the barotropic stepping
539 !! [H ~> m or kg m-2].
540 real, dimension(:,:), pointer :: taux_bot !< The zonal bottom frictional stress from
541 !! ocean to the seafloor [R L Z T-2 ~> Pa].
542 real, dimension(:,:), pointer :: tauy_bot !< The meridional bottom frictional stress
543 !! from ocean to the seafloor [R L Z T-2 ~> Pa].
544 real, dimension(:,:,:), pointer :: uh0 !< The zonal layer transports at reference
545 !! velocities [H L2 T-1 ~> m3 s-1 or kg s-1].
546 real, dimension(:,:,:), pointer :: u_uh0 !< The velocities used to calculate
547 !! uh0 [L T-1 ~> m s-1]
548 real, dimension(:,:,:), pointer :: vh0 !< The zonal layer transports at reference
549 !! velocities [H L2 T-1 ~> m3 s-1 or kg s-1].
550 real, dimension(:,:,:), pointer :: v_vh0 !< The velocities used to calculate
551 !! vh0 [L T-1 ~> m s-1]
552 real, dimension(SZI_(G),SZJ_(G)), optional, intent(out) :: etaav !< The free surface height or column mass
553 !! averaged over the barotropic integration [H ~> m or kg m-2].
554
555 ! Local variables
556 real :: ubt_cor(szib_(g),szj_(g)) ! The barotropic velocities that had been
557 real :: vbt_cor(szi_(g),szjb_(g)) ! used to calculate the input Coriolis
558 ! terms [L T-1 ~> m s-1].
559 real :: wt_u(szib_(g),szj_(g),szk_(gv)) ! wt_u and wt_v are the
560 real :: wt_v(szi_(g),szjb_(g),szk_(gv)) ! normalized weights to
561 ! be used in calculating barotropic velocities, possibly with
562 ! sums less than one due to viscous losses [nondim]
563 real :: iwt_u_tot(szib_(g),szj_(g)) ! Iwt_u_tot and Iwt_v_tot are the
564 real :: iwt_v_tot(szi_(g),szjb_(g)) ! inverses of wt_u and wt_v vertical integrals,
565 ! used to normalize wt_u and wt_v [nondim]
566 real, dimension(SZIB_(G),SZJ_(G)) :: &
567 av_rem_u, & ! The weighted average of visc_rem_u [nondim]
568 tmp_u, & ! A temporary array at u points [L T-2 ~> m s-2] or [nondim]
569 ubt_st, & ! The zonal barotropic velocity at the start of timestep [L T-1 ~> m s-1].
570 ubt_wtd, & ! A weighted sum used to find the filtered final ubt [L T-1 ~> m s-1].
571 pfu_avg, & ! The average zonal barotropic pressure gradient force [L T-2 ~> m s-2].
572 coru_avg, & ! The average zonal barotropic Coriolis acceleration [L T-2 ~> m s-2].
573 ldu_avg, & ! The average zonal barotropic linear wave drag acceleration [L T-2 ~> m s-2].
574 ubt_dt ! The zonal barotropic velocity tendency [L T-2 ~> m s-2].
575 real, dimension(SZI_(G),SZJB_(G)) :: &
576 av_rem_v, & ! The weighted average of visc_rem_v [nondim]
577 tmp_v, & ! A temporary array at v points [L T-2 ~> m s-2] or [nondim]
578 vbt_st, & ! The meridional barotropic velocity at the start of timestep [L T-1 ~> m s-1].
579 vbt_wtd, & ! A weighted sum used to find the filtered final vbt [L T-1 ~> m s-1].
580 pfv_avg, & ! The average meridional barotropic pressure gradient force [L T-2 ~> m s-2].
581 corv_avg, & ! The average meridional barotropic Coriolis acceleration [L T-2 ~> m s-2].
582 ldv_avg, & ! The average meridional barotropic linear wave drag acceleration [L T-2 ~> m s-2].
583 vbt_dt ! The meridional barotropic velocity tendency [L T-2 ~> m s-2].
584 real, dimension(SZI_(G),SZJ_(G)) :: &
585 tmp_h, & ! A temporary array at h points [nondim]
586 e_anom ! The anomaly in the sea surface height or column mass
587 ! averaged between the beginning and end of the time step,
588 ! relative to eta_PF, with SAL effects included [H ~> m or kg m-2].
589
590 ! These are always allocated with symmetric memory and wide halos.
591 real :: q(szibw_(cs),szjbw_(cs)) ! A pseudo potential vorticity [T-1 H-1 ~> s-1 m-1 or m2 s-1 kg-1]
592 real, dimension(SZIBW_(CS),SZJW_(CS)) :: &
593 ubt, & ! The zonal barotropic velocity [L T-1 ~> m s-1].
594 bt_rem_u, & ! The fraction of the barotropic zonal velocity that remains
595 ! after a time step, the remainder being lost to bottom drag [nondim].
596 ! bt_rem_u is between 0 and 1.
597 bt_force_u, & ! The vertical average of all of the u-accelerations that are
598 ! not explicitly included in the barotropic equation [L T-2 ~> m s-2].
599 u_accel_bt, & ! The difference between the zonal acceleration from the
600 ! barotropic calculation and BT_force_u [L T-2 ~> m s-2].
601 uhbt, & ! The zonal barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].
602 uhbt0, & ! The difference between the sum of the layer zonal thickness
603 ! fluxes and the barotropic thickness flux using the same
604 ! velocity [H L2 T-1 ~> m3 s-1 or kg s-1].
605 cor_ref_u, & ! The zonal barotropic Coriolis acceleration due
606 ! to the reference velocities [L T-2 ~> m s-2].
607 rayleigh_u, & ! A Rayleigh drag timescale operating at u-points for drag parameterizations
608 ! that introduced directly into the barotropic solver rather than coming in via
609 ! the visc_rem_u arrays from the layered equations [T-1 ~> s-1].
610 ! This is nonzero mostly for a barotropic tidal body drag.
611 dcor_u, & ! An averaged total thickness at u points [H ~> m or kg m-2].
612 datu ! Basin depth at u-velocity grid points times the y-grid
613 ! spacing [H L ~> m2 or kg m-1].
614 real, dimension(SZIW_(CS),SZJBW_(CS)) :: &
615 vbt, & ! The meridional barotropic velocity [L T-1 ~> m s-1].
616 bt_rem_v, & ! The fraction of the barotropic meridional velocity that
617 ! remains after a time step, the rest being lost to bottom
618 ! drag [nondim]. bt_rem_v is between 0 and 1.
619 bt_force_v, & ! The vertical average of all of the v-accelerations that are
620 ! not explicitly included in the barotropic equation [L T-2 ~> m s-2].
621 v_accel_bt, & ! The difference between the meridional acceleration from the
622 ! barotropic calculation and BT_force_v [L T-2 ~> m s-2].
623 vhbt, & ! The meridional barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].
624 vhbt0, & ! The difference between the sum of the layer meridional
625 ! thickness fluxes and the barotropic thickness flux using
626 ! the same velocities [H L2 T-1 ~> m3 s-1 or kg s-1].
627 cor_ref_v, & ! The meridional barotropic Coriolis acceleration due
628 ! to the reference velocities [L T-2 ~> m s-2].
629 rayleigh_v, & ! A Rayleigh drag timescale operating at v-points for drag parameterizations
630 ! that introduced directly into the barotropic solver rather than coming
631 ! in via the visc_rem_v arrays from the layered equations [T-1 ~> s-1].
632 ! This is nonzero mostly for a barotropic tidal body drag.
633 dcor_v, & ! An averaged total thickness at v points [H ~> m or kg m-2].
634 datv ! Basin depth at v-velocity grid points times the x-grid
635 ! spacing [H L ~> m2 or kg m-1].
636 real, dimension(4,SZIBW_(CS),SZJW_(CS)) :: &
637 f_4_u !< The terms giving the contribution to the Coriolis acceleration at a zonal
638 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
639 !! These are the products of thicknesses at v points and appropriately staggered
640 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
641 !! they are approximately f / 4. The 4 values on the innermost loop are for
642 !! v-velocities to the southwest, southeast, northwest and northeast.
643 real, dimension(4,SZIW_(CS),SZJBW_(CS)) :: &
644 f_4_v !< The terms giving the contribution to the Coriolis acceleration at a meridional
645 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
646 !! These are the products of thicknesses at u points and appropriately staggered
647 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
648 !! they are approximately f / 4. The 4 values on the innermost loop are for
649 !! u-velocities to the southwest, southeast, northwest and northeast.
650 real, dimension(:,:,:), pointer :: ufilt, vfilt
651 ! Filtered velocities from the output of streaming filters [L T-1 ~> m s-1]
652 real, dimension(SZIB_(G),SZJ_(G)) :: drag_u
653 ! The zonal acceleration due to frequency-dependent drag [L T-2 ~> m s-2]
654 real, dimension(SZI_(G),SZJB_(G)) :: drag_v
655 ! The meridional acceleration due to frequency-dependent drag [L T-2 ~> m s-2]
656 real, target, dimension(SZIW_(CS),SZJW_(CS)) :: &
657 eta ! The barotropic free surface height anomaly or column mass
658 ! anomaly [H ~> m or kg m-2]
659 real, dimension(SZIW_(CS),SZJW_(CS)) :: &
660 eta_sum, & ! eta summed across the timesteps [H ~> m or kg m-2].
661 eta_wtd, & ! A weighted estimate used to calculate eta_out [H ~> m or kg m-2].
662 eta_ic, & ! A local copy of the initial 2-D eta field (eta_in) [H ~> m or kg m-2]
663 eta_pf, & ! A local copy of the 2-D eta field (either SSH anomaly or
664 ! column mass anomaly) that was used to calculate the input
665 ! pressure gradient accelerations [H ~> m or kg m-2].
666 eta_pf_1, & ! The initial value of eta_PF, when interp_eta_PF is
667 ! true [H ~> m or kg m-2].
668 d_eta_pf, & ! The change in eta_PF over the barotropic time stepping when
669 ! interp_eta_PF is true [H ~> m or kg m-2].
670 gtot_e, & ! gtot_X is the effective total reduced gravity used to relate
671 gtot_w, & ! free surface height deviations to pressure forces (including
672 gtot_n, & ! GFS and baroclinic contributions) in the barotropic momentum
673 gtot_s, & ! equations half a grid-point in the X-direction (X is N, S, E, or W)
674 ! from the thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
675 ! (See Hallberg, J Comp Phys 1997 for a discussion.)
676 eta_src, & ! The source of eta per barotropic timestep [H ~> m or kg m-2].
677 spv_col_avg, & ! The column average specific volume [R-1 ~> m3 kg-1]
678 dyn_coef_eta ! The coefficient relating the changes in eta to the
679 ! dynamic surface pressure under rigid ice
680 ! [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
681 type(local_bt_cont_u_type), dimension(SZIBW_(CS),SZJW_(CS)) :: &
682 btcl_u ! A repackaged version of the u-point information in BT_cont.
683 type(local_bt_cont_v_type), dimension(SZIW_(CS),SZJBW_(CS)) :: &
684 btcl_v ! A repackaged version of the v-point information in BT_cont.
685 ! End of wide-sized variables.
686
687 real :: visc_rem ! A work variable that may equal visc_rem_[uv] [nondim]
688 real :: dtbt ! The barotropic time step [T ~> s].
689 real :: idt ! The inverse of dt [T-1 ~> s-1].
690 real :: det_de ! The partial derivative due to self-attraction and loading
691 ! of the reference geopotential with the sea surface height [nondim].
692 ! This is typically ~0.09 or less.
693 real :: dgeo_de ! The constant of proportionality between geopotential and sea surface height
694 ! [nondim]. It is of order 1, but for stability this may be made larger than
695 ! the physical problem would suggest.
696 real :: dgeo_de_obc ! The value of dgeo_de to be used with Flather open boundary conditions [nondim].
697 real :: instep ! The inverse of the number of barotropic time steps to take [nondim].
698 integer :: nstep ! The number of barotropic time steps to take.
699 real :: htot_avg ! The average total thickness of the tracer columns adjacent to a
700 ! velocity point [H ~> m or kg m-2]
701 logical :: use_bt_cont, find_etaav
702 logical :: integral_bt_cont ! If true, update the barotropic continuity equation directly
703 ! from the initial condition using the time-integrated barotropic velocity.
704 logical :: ice_is_rigid, nonblock_setup, interp_eta_pf
705 logical :: add_uh0
706
707 real :: dyn_coef_max ! The maximum stable value of dyn_coef_eta
708 ! [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
709 real :: ice_strength = 0.0 ! The effective strength of the ice [L2 Z-1 T-2 ~> m s-2].
710 real :: h_to_z ! A local unit conversion factor used with rigid ice [Z H-1 ~> nondim or m3 kg-1]
711 real :: idt_max2 ! The squared inverse of the local maximum stable
712 ! barotropic time step [T-2 ~> s-2].
713 real :: h_min_dyn ! The minimum depth to use in limiting the size of the
714 ! dynamic surface pressure for stability [H ~> m or kg m-2].
715 real :: h_eff_dx2 ! The effective total thickness divided by the grid spacing
716 ! squared [H L-2 ~> m-1 or kg m-4].
717 real :: u_max_cor, v_max_cor ! The maximum corrective velocities [L T-1 ~> m s-1].
718 real :: uint_cor, vint_cor ! The maximum time-integrated corrective velocities [L ~> m].
719 real :: htot ! The total thickness [H ~> m or kg m-2].
720 real :: eta_cor_max ! The maximum fluid that can be added as a correction to eta [H ~> m or kg m-2].
721 real :: accel_underflow ! An acceleration that is so small it should be zeroed out [L T-2 ~> m s-2].
722 real :: h_a_neglect ! A cell volume or mass that is so small it is usually lost
723 ! in roundoff and can be neglected [H L2 ~> m3 or kg].
724
725 real, allocatable :: wt_vel(:) ! The raw or relative weights of each of the barotropic timesteps
726 ! in determining the average velocities [nondim]
727 real, allocatable :: wt_eta(:) ! The raw or relative weights of each of the barotropic timesteps
728 ! in determining the average eta [nondim]
729 real, allocatable :: wt_accel(:) ! The raw or relative weights of each of the barotropic timesteps
730 ! in determining the average accelerations [nondim]
731 real, allocatable :: wt_trans(:) ! The raw or relative weights of each of the barotropic timesteps
732 ! in determining the average transports [nondim]
733 real, allocatable :: wt_accel2(:) ! A potentially un-normalized copy of wt_accel [nondim]
734 real :: sum_wt_vel ! The sum of the raw weights used to find average velocities [nondim]
735 real :: sum_wt_eta ! The sum of the raw weights used to find average eta [nondim]
736 real :: sum_wt_accel ! The sum of the raw weights used to find average accelerations [nondim]
737 real :: sum_wt_trans ! The sum of the raw weights used to find average transports [nondim]
738 real :: i_sum_wt_vel ! The inverse of the sum of the raw weights used to find average velocities [nondim]
739 real :: i_sum_wt_eta ! The inverse of the sum of the raw weights used to find eta [nondim]
740 real :: i_sum_wt_accel ! The inverse of the sum of the raw weights used to find average accelerations [nondim]
741 real :: i_sum_wt_trans ! The inverse of the sum of the raw weights used to find average transports [nondim]
742 real :: dt_filt ! The half-width of the barotropic filter [T ~> s].
743 integer :: nfilter
744
745 logical :: apply_obcs, apply_obc_flather
746 type(memory_size_type) :: ms
747 character(len=200) :: mesg
748 integer :: stencil ! The stencil size of the algorithm, often 1 or 2.
749 integer :: isvf, ievf, jsvf, jevf, num_cycles
750 integer :: i, j, k, n
751 integer :: is, ie, js, je, nz, isq, ieq, jsq, jeq
752 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
753
754 if (.not.cs%module_is_initialized) call mom_error(fatal, &
755 "btstep: Module MOM_barotropic must be initialized before it is used.")
756
757 if (.not.cs%split) return
758 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
759 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
760 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
761 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
762 ms%isdw = cs%isdw ; ms%iedw = cs%iedw ; ms%jsdw = cs%jsdw ; ms%jedw = cs%jedw
763 h_a_neglect = gv%H_subroundoff * (1.0 * us%m_to_L**2)
764
765 idt = 1.0 / dt
766 accel_underflow = cs%vel_underflow * idt
767
768 use_bt_cont = associated(bt_cont)
769 integral_bt_cont = use_bt_cont .and. cs%integral_BT_cont
770
771 interp_eta_pf = associated(eta_pf_start)
772
773 ! Figure out the fullest arrays that could be updated.
774 stencil = max(1, cs%min_stencil)
775 if ((.not.use_bt_cont) .and. cs%Nonlinear_continuity .and. &
776 (cs%Nonlin_cont_update_period > 0)) stencil = max(2, cs%min_stencil)
777
778 find_etaav = present(etaav)
779
780 add_uh0 = associated(uh0)
781 if (add_uh0 .and. .not.(associated(vh0) .and. associated(u_uh0) .and. &
782 associated(v_vh0))) call mom_error(fatal, &
783 "btstep: vh0, u_uh0, and v_vh0 must be associated if uh0 is used.")
784
785 ! This can be changed to try to optimize the performance.
786 nonblock_setup = g%nonblocking_updates
787
788 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
789
790 apply_obc_flather = .false.
791 apply_obcs = .false.
792 if (associated(obc)) then
793 apply_obc_flather = open_boundary_query(obc, apply_flather_obc=.true.)
794 apply_obcs = open_boundary_query(obc, apply_specified_obc=.true.) .or. &
795 apply_obc_flather .or. open_boundary_query(obc, apply_open_obc=.true.)
796 endif
797
798 num_cycles = 1
799 if (cs%use_wide_halos) &
800 num_cycles = min((is-cs%isdw) / stencil, (js-cs%jsdw) / stencil)
801 isvf = is - (num_cycles-1)*stencil ; ievf = ie + (num_cycles-1)*stencil
802 jsvf = js - (num_cycles-1)*stencil ; jevf = je + (num_cycles-1)*stencil
803
804 nstep = ceiling(dt/cs%dtbt - 0.0001)
805 if (is_root_pe() .and. ((nstep /= cs%nstep_last) .or. cs%debug)) then
806 write(mesg,'("btstep is using a dynamic barotropic timestep of ", ES12.6, &
807 & " seconds, max ", ES12.6, ".")') (us%T_to_s*dt/nstep), us%T_to_s*cs%dtbt_max
808 call mom_mesg(mesg, 3)
809 endif
810 cs%nstep_last = nstep
811
812 ! Set the actual barotropic time step.
813 instep = 1.0 / real(nstep)
814 dtbt = dt * instep
815
816 !--- begin setup for group halo update
817 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
818 if (.not. cs%linearized_BT_PV) then
819 call create_group_pass(cs%pass_q_DCor, q, cs%BT_Domain, to_all, position=corner)
820 call create_group_pass(cs%pass_q_DCor, dcor_u, dcor_v, cs%BT_Domain, &
821 to_all+scalar_pair)
822 endif
823 if ((isq > is-1) .or. (jsq > js-1)) &
824 call create_group_pass(cs%pass_tmp_uv, tmp_u, tmp_v, g%Domain)
825 call create_group_pass(cs%pass_gtot, gtot_e, gtot_n, cs%BT_Domain, &
826 to_all+scalar_pair, agrid)
827 call create_group_pass(cs%pass_gtot, gtot_w, gtot_s, cs%BT_Domain, &
828 to_all+scalar_pair, agrid)
829
830 if (cs%dynamic_psurf) &
831 call create_group_pass(cs%pass_eta_bt_rem, dyn_coef_eta, cs%BT_Domain)
832 if (interp_eta_pf) then
833 call create_group_pass(cs%pass_eta_bt_rem, eta_pf_1, cs%BT_Domain)
834 call create_group_pass(cs%pass_eta_bt_rem, d_eta_pf, cs%BT_Domain)
835 else
836 call create_group_pass(cs%pass_eta_bt_rem, eta_pf, cs%BT_Domain)
837 endif
838 if (integral_bt_cont) &
839 call create_group_pass(cs%pass_eta_bt_rem, eta_ic, cs%BT_Domain)
840 call create_group_pass(cs%pass_eta_bt_rem, eta_src, cs%BT_Domain)
841
842 call create_group_pass(cs%pass_eta_bt_rem, bt_rem_u, bt_rem_v, &
843 cs%BT_Domain, to_all+scalar_pair)
844 if (cs%linear_wave_drag) &
845 call create_group_pass(cs%pass_eta_bt_rem, rayleigh_u, rayleigh_v, &
846 cs%BT_Domain, to_all+scalar_pair)
847
848 ! The following halo update is not needed without wide halos. RWH
849 if (((g%isd > cs%isdw) .or. (g%jsd > cs%jsdw)) .or. (isq <= is-1) .or. (jsq <= js-1)) &
850 call create_group_pass(cs%pass_force_hbt0_Cor_ref, bt_force_u, bt_force_v, cs%BT_Domain)
851 if (add_uh0) call create_group_pass(cs%pass_force_hbt0_Cor_ref, uhbt0, vhbt0, cs%BT_Domain)
852 call create_group_pass(cs%pass_force_hbt0_Cor_ref, cor_ref_u, cor_ref_v, cs%BT_Domain)
853 if (.not. use_bt_cont) then
854 call create_group_pass(cs%pass_Dat_uv, datu, datv, cs%BT_Domain, to_all+scalar_pair)
855 endif
856 if (apply_obc_flather .and. .not.gv%Boussinesq) &
857 call create_group_pass(cs%pass_SpV_avg, spv_col_avg, cs%BT_domain)
858
859 call create_group_pass(cs%pass_ubt_Cor, ubt_cor, vbt_cor, g%Domain)
860 ! These passes occur at the end of the routine, as data is being readied to
861 ! share with the main part of the MOM6 code.
862 if (find_etaav) then
863 call create_group_pass(cs%pass_etaav, etaav, g%Domain)
864 endif
865 call create_group_pass(cs%pass_e_anom, e_anom, g%Domain)
866 call create_group_pass(cs%pass_ubta_uhbta, cs%ubtav, cs%vbtav, g%Domain)
867 call create_group_pass(cs%pass_ubta_uhbta, uhbtav, vhbtav, g%Domain)
868
869 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
870!--- end setup for group halo update
871
872! Calculate the constant coefficients for the Coriolis force terms in the
873! barotropic momentum equations. This has to be done quite early to start
874! the halo update that needs to be completed before the next calculations.
875 if (cs%linearized_BT_PV) then
876 !$OMP parallel do default(shared)
877 do j=jsvf-2,jevf+1 ; do i=isvf-2,ievf+1
878 q(i,j) = cs%q_D(i,j)
879 enddo ; enddo
880 !$OMP parallel do default(shared)
881 do j=jsvf-1,jevf+1 ; do i=isvf-2,ievf+1
882 dcor_u(i,j) = cs%D_u_Cor(i,j)
883 enddo ; enddo
884 !$OMP parallel do default(shared)
885 do j=jsvf-2,jevf+1 ; do i=isvf-1,ievf+1
886 dcor_v(i,j) = cs%D_v_Cor(i,j)
887 enddo ; enddo
888 else
889 q(:,:) = 0.0 ; dcor_u(:,:) = 0.0 ; dcor_v(:,:) = 0.0
890 if (gv%Boussinesq) then
891 !$OMP parallel do default(shared)
892 do j=js,je ; do i=is-1,ie
893 dcor_u(i,j) = 0.5 * (max(gv%Z_to_H*g%bathyT(i+1,j) + eta_in(i+1,j), 0.0) + &
894 max(gv%Z_to_H*g%bathyT(i,j) + eta_in(i,j), 0.0) )
895 enddo ; enddo
896 if (cs%interior_OBC_PV .and. cs%BT_OBC%u_OBCs_on_PE) then
897 !$OMP parallel do default(shared)
898 do j = max(js,cs%BT_OBC%js_u_W_obc), min(je,cs%BT_OBC%je_u_W_obc)
899 do i = max(is-1,cs%BT_OBC%Is_u_W_obc), min(ie,cs%BT_OBC%Ie_u_W_obc)
900 if (cs%BT_OBC%u_OBC_type(i,j) < 0) & ! Western boundary condition
901 dcor_u(i,j) = max(gv%Z_to_H*g%bathyT(i+1,j) + eta_in(i+1,j), 0.0)
902 enddo
903 enddo
904 !$OMP parallel do default(shared)
905 do j = max(js,cs%BT_OBC%js_u_E_obc), min(je,cs%BT_OBC%je_u_E_obc)
906 do i = max(is-1,cs%BT_OBC%Is_u_E_obc), min(ie,cs%BT_OBC%Ie_u_E_obc)
907 if (cs%BT_OBC%u_OBC_type(i,j) > 0) & ! Eastern boundary condition
908 dcor_u(i,j) = max(gv%Z_to_H*g%bathyT(i,j) + eta_in(i,j), 0.0)
909 enddo
910 enddo
911 endif
912
913 !$OMP parallel do default(shared)
914 do j=js-1,je ; do i=is,ie
915 dcor_v(i,j) = 0.5 * (max(gv%Z_to_H*g%bathyT(i,j+1) + eta_in(i+1,j), 0.0) + &
916 max(gv%Z_to_H*g%bathyT(i,j) + eta_in(i,j), 0.0) )
917 enddo ; enddo
918 if (cs%interior_OBC_PV .and. cs%BT_OBC%v_OBCs_on_PE) then
919 !$OMP parallel do default(shared)
920 do j = max(js,cs%BT_OBC%js_v_S_obc), min(je,cs%BT_OBC%je_v_S_obc)
921 do i = max(is-1,cs%BT_OBC%Is_v_S_obc), min(ie,cs%BT_OBC%Ie_v_S_obc)
922 if (cs%BT_OBC%v_OBC_type(i,j) < 0) & ! Southern boundary condition
923 dcor_v(i,j) = max(gv%Z_to_H*g%bathyT(i,j+1) + eta_in(i,j+1), 0.0)
924 enddo
925 enddo
926 !$OMP parallel do default(shared)
927 do j = max(js,cs%BT_OBC%js_v_N_obc), min(je,cs%BT_OBC%je_v_N_obc)
928 do i = max(is-1,cs%BT_OBC%Is_v_N_obc), min(ie,cs%BT_OBC%Ie_v_N_obc)
929 if (cs%BT_OBC%v_OBC_type(i,j) > 0) & ! Northern boundary condition
930 dcor_v(i,j) = max(gv%Z_to_H*g%bathyT(i,j) + eta_in(i,j), 0.0)
931 enddo
932 enddo
933 endif
934 !$OMP parallel do default(shared)
935 do j=js-1,je ; do i=is-1,ie
936 q(i,j) = 0.25 * (cs%BT_Coriolis_scale * g%CoriolisBu(i,j)) * &
937 ((cs%q_wt(1,i,j) + cs%q_wt(4,i,j)) + (cs%q_wt(2,i,j) + cs%q_wt(3,i,j))) / &
938 (max(((cs%q_wt(1,i,j) * max(gv%Z_to_H*g%bathyT(i,j) + eta_in(i,j), 0.0)) + &
939 (cs%q_wt(4,i,j) * max(gv%Z_to_H*g%bathyT(i+1,j+1) + eta_in(i+1,j+1), 0.0))) + &
940 ((cs%q_wt(2,i,j) * max(gv%Z_to_H*g%bathyT(i+1,j) + eta_in(i+1,j), 0.0)) + &
941 (cs%q_wt(3,i,j) * max(gv%Z_to_H*g%bathyT(i,j+1) + eta_in(i,j+1), 0.0))), h_a_neglect) )
942 enddo ; enddo
943 else ! Non-Boussinesq
944 !$OMP parallel do default(shared)
945 do j=js,je ; do i=is-1,ie
946 dcor_u(i,j) = 0.5 * (eta_in(i+1,j) + eta_in(i,j))
947 enddo ; enddo
948 if (cs%interior_OBC_PV .and. cs%BT_OBC%u_OBCs_on_PE) then
949 !$OMP parallel do default(shared)
950 do j = max(js,cs%BT_OBC%js_u_W_obc), min(je,cs%BT_OBC%je_u_W_obc)
951 do i = max(is-1,cs%BT_OBC%Is_u_W_obc), min(ie,cs%BT_OBC%Ie_u_W_obc)
952 if (cs%BT_OBC%u_OBC_type(i,j) < 0) dcor_u(i,j) = eta_in(i+1,j) ! Western boundary condition
953 enddo
954 enddo
955 !$OMP parallel do default(shared)
956 do j = max(js,cs%BT_OBC%js_u_E_obc), min(je,cs%BT_OBC%je_u_E_obc)
957 do i = max(is-1,cs%BT_OBC%Is_u_E_obc), min(ie,cs%BT_OBC%Ie_u_E_obc)
958 if (cs%BT_OBC%u_OBC_type(i,j) > 0) dcor_u(i,j) = eta_in(i,j) ! Eastern boundary condition
959 enddo
960 enddo
961 endif
962
963 !$OMP parallel do default(shared)
964 do j=js-1,je ; do i=is,ie
965 dcor_v(i,j) = 0.5 * (eta_in(i,j+1) + eta_in(i,j))
966 enddo ; enddo
967 if (cs%interior_OBC_PV .and. cs%BT_OBC%v_OBCs_on_PE) then
968 !$OMP parallel do default(shared)
969 do j = max(js,cs%BT_OBC%js_v_S_obc), min(je,cs%BT_OBC%je_v_S_obc)
970 do i = max(is-1,cs%BT_OBC%Is_v_S_obc), min(ie,cs%BT_OBC%Ie_v_S_obc)
971 if (cs%BT_OBC%v_OBC_type(i,j) < 0) dcor_v(i,j) = eta_in(i,j+1) ! Southern boundary condition
972 enddo
973 enddo
974 !$OMP parallel do default(shared)
975 do j = max(js,cs%BT_OBC%js_v_N_obc), min(je,cs%BT_OBC%je_v_N_obc)
976 do i = max(is-1,cs%BT_OBC%Is_v_N_obc), min(ie,cs%BT_OBC%Ie_v_N_obc)
977 if (cs%BT_OBC%v_OBC_type(i,j) > 0) dcor_v(i,j) = eta_in(i,j) ! Northern boundary condition
978 enddo
979 enddo
980 endif
981
982 !$OMP parallel do default(shared)
983 do j=js-1,je ; do i=is-1,ie
984 q(i,j) = 0.25 * (cs%BT_Coriolis_scale * g%CoriolisBu(i,j)) * &
985 ((cs%q_wt(1,i,j) + cs%q_wt(4,i,j)) + (cs%q_wt(2,i,j) + cs%q_wt(3,i,j))) / &
986 (max(((cs%q_wt(1,i,j) * eta_in(i,j)) + (cs%q_wt(4,i,j) * eta_in(i+1,j+1))) + &
987 ((cs%q_wt(2,i,j) * eta_in(i+1,j)) + (cs%q_wt(3,i,j) * eta_in(i,j+1))), h_a_neglect) )
988 enddo ; enddo
989 endif
990
991 ! With very wide halos, q and D need to be calculated on the available data
992 ! domain and then updated onto the full computational domain.
993 ! These calculations can be done almost immediately, but the halo updates
994 ! must be done before the [abcd]mer and [abcd]zon are calculated.
995 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
996 if (nonblock_setup) then
997 call start_group_pass(cs%pass_q_DCor, cs%BT_Domain, clock=id_clock_pass_pre)
998 else
999 call do_group_pass(cs%pass_q_DCor, cs%BT_Domain, clock=id_clock_pass_pre)
1000 endif
1001 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1002 endif
1003
1004 ! Zero out various wide-halo arrays.
1005 !$OMP parallel do default(shared)
1006 do j=cs%jsdw,cs%jedw ; do i=cs%isdw,cs%iedw
1007 gtot_e(i,j) = 0.0 ; gtot_w(i,j) = 0.0
1008 gtot_n(i,j) = 0.0 ; gtot_s(i,j) = 0.0
1009 eta(i,j) = 0.0
1010 eta_pf(i,j) = 0.0
1011 if (interp_eta_pf) then
1012 eta_pf_1(i,j) = 0.0 ; d_eta_pf(i,j) = 0.0
1013 endif
1014 if (integral_bt_cont) then
1015 eta_ic(i,j) = 0.0
1016 endif
1017 if (cs%dynamic_psurf) dyn_coef_eta(i,j) = 0.0
1018 enddo ; enddo
1019 ! The halo regions of various arrays need to be initialized to
1020 ! non-NaNs in case the neighboring domains are not part of the ocean.
1021 ! Otherwise a halo update later on fills in the correct values.
1022 !$OMP parallel do default(shared)
1023 do j=cs%jsdw,cs%jedw ; do i=cs%isdw-1,cs%iedw
1024 cor_ref_u(i,j) = 0.0 ; bt_force_u(i,j) = 0.0 ; ubt(i,j) = 0.0
1025 datu(i,j) = 0.0 ; bt_rem_u(i,j) = 0.0 ; uhbt0(i,j) = 0.0
1026 enddo ; enddo
1027 !$OMP parallel do default(shared)
1028 do j=cs%jsdw-1,cs%jedw ; do i=cs%isdw,cs%iedw
1029 cor_ref_v(i,j) = 0.0 ; bt_force_v(i,j) = 0.0 ; vbt(i,j) = 0.0
1030 datv(i,j) = 0.0 ; bt_rem_v(i,j) = 0.0 ; vhbt0(i,j) = 0.0
1031 enddo ; enddo
1032
1033 if (apply_obcs) then
1034 spv_col_avg(:,:) = 0.0
1035 if (apply_obc_flather .and. .not.gv%Boussinesq) then
1036 ! Copy the column average specific volumes into a wide halo array
1037 !$OMP parallel do default(shared)
1038 do j=js,je ; do i=is,ie
1039 spv_col_avg(i,j) = spv_avg(i,j)
1040 enddo ; enddo
1041 if (nonblock_setup) then
1042 call start_group_pass(cs%pass_SpV_avg, cs%BT_domain)
1043 else
1044 call do_group_pass(cs%pass_SpV_avg, cs%BT_domain)
1045 endif
1046 endif
1047 endif
1048
1049 if (cs%linear_wave_drag) then
1050 !$OMP parallel do default(shared)
1051 do j=cs%jsdw,cs%jedw ; do i=cs%isdw-1,cs%iedw
1052 rayleigh_u(i,j) = 0.0
1053 enddo ; enddo
1054 !$OMP parallel do default(shared)
1055 do j=cs%jsdw-1,cs%jedw ; do i=cs%isdw,cs%iedw
1056 rayleigh_v(i,j) = 0.0
1057 enddo ; enddo
1058 endif
1059
1060 ! Copy input arrays into their wide-halo counterparts.
1061 if (interp_eta_pf) then
1062 !$OMP parallel do default(shared)
1063 do j=g%jsd,g%jed ; do i=g%isd,g%ied ! Was "do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1" but doing so breaks OBC. Not sure why?
1064 eta(i,j) = eta_in(i,j)
1065 eta_pf_1(i,j) = eta_pf_start(i,j)
1066 d_eta_pf(i,j) = eta_pf_in(i,j) - eta_pf_start(i,j)
1067 enddo ; enddo
1068 else
1069 !$OMP parallel do default(shared)
1070 do j=g%jsd,g%jed ; do i=g%isd,g%ied !: Was "do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1" but doing so breaks OBC. Not sure why?
1071 eta(i,j) = eta_in(i,j)
1072 eta_pf(i,j) = eta_pf_in(i,j)
1073 enddo ; enddo
1074 endif
1075 if (integral_bt_cont) then
1076 !$OMP parallel do default(shared)
1077 do j=g%jsd,g%jed ; do i=g%isd,g%ied
1078 eta_ic(i,j) = eta_in(i,j)
1079 enddo ; enddo
1080 endif
1081
1082 !$OMP parallel do default(shared) private(visc_rem)
1083 do k=1,nz ; do j=js,je ; do i=is-1,ie
1084 ! rem needs to be greater than visc_rem_u and 1-Instep/visc_rem_u.
1085 ! The 0.5 below is just for safety.
1086 ! NOTE: subroundoff is a negligible value used to prevent division by zero.
1087 ! When 1-0.5*Instep/visc_rem exceeds visc_rem, the subroundoff is too small
1088 ! to modify the significand. When visc_rem is small, the max() operators
1089 ! select visc_rem or 0. So subroundoff cannot impact the final value.
1090 visc_rem = min(visc_rem_u(i,j,k), 1.)
1091 visc_rem = max(visc_rem, 1. - 0.5 * instep / (visc_rem + subroundoff))
1092 visc_rem = max(visc_rem, 0.)
1093 wt_u(i,j,k) = cs%frhatu(i,j,k) * visc_rem
1094 enddo ; enddo ; enddo
1095 !$OMP parallel do default(shared) private(visc_rem)
1096 do k=1,nz ; do j=js-1,je ; do i=is,ie
1097 ! As above, rem must be greater than visc_rem_v and 1-Instep/visc_rem_v.
1098 visc_rem = min(visc_rem_v(i,j,k), 1.)
1099 visc_rem = max(visc_rem, 1. - 0.5 * instep / (visc_rem + subroundoff))
1100 visc_rem = max(visc_rem, 0.)
1101 wt_v(i,j,k) = cs%frhatv(i,j,k) * visc_rem
1102 enddo ; enddo ; enddo
1103
1104 if (.not. cs%wt_uv_bug) then
1105 do j=js,je ; do i=is-1,ie ; iwt_u_tot(i,j) = wt_u(i,j,1) ; enddo ; enddo
1106 do k=2,nz ; do j=js,je ; do i=is-1,ie
1107 iwt_u_tot(i,j) = iwt_u_tot(i,j) + wt_u(i,j,k)
1108 enddo ; enddo ; enddo
1109 do j=js,je ; do i=is-1,ie
1110 if (abs(iwt_u_tot(i,j)) > 0.0 ) iwt_u_tot(i,j) = g%mask2dCu(i,j) / iwt_u_tot(i,j)
1111 enddo ; enddo
1112 do k=1,nz ; do j=js,je ; do i=is-1,ie
1113 wt_u(i,j,k) = wt_u(i,j,k) * iwt_u_tot(i,j)
1114 enddo ; enddo ; enddo
1115
1116 do j=js-1,je ; do i=is,ie ; iwt_v_tot(i,j) = wt_v(i,j,1) ; enddo ; enddo
1117 do k=2,nz ; do j=js-1,je ; do i=is,ie
1118 iwt_v_tot(i,j) = iwt_v_tot(i,j) + wt_v(i,j,k)
1119 enddo ; enddo ; enddo
1120 do j=js-1,je ; do i=is,ie
1121 if (abs(iwt_v_tot(i,j)) > 0.0 ) iwt_v_tot(i,j) = g%mask2dCv(i,j) / iwt_v_tot(i,j)
1122 enddo ; enddo
1123 do k=1,nz ; do j=js-1,je ; do i=is,ie
1124 wt_v(i,j,k) = wt_v(i,j,k) * iwt_v_tot(i,j)
1125 enddo ; enddo ; enddo
1126 endif
1127
1128 ! Use u_Cor and v_Cor as the reference values for the Coriolis terms,
1129 ! including the viscous remnant.
1130 !$OMP parallel do default(shared)
1131 do j=js-1,je+1 ; do i=is-1,ie ; ubt_cor(i,j) = 0.0 ; enddo ; enddo
1132 !$OMP parallel do default(shared)
1133 do j=js-1,je ; do i=is-1,ie+1 ; vbt_cor(i,j) = 0.0 ; enddo ; enddo
1134 !$OMP parallel do default(shared)
1135 do j=js,je ; do k=1,nz ; do i=is-1,ie
1136 ubt_cor(i,j) = ubt_cor(i,j) + wt_u(i,j,k) * u_cor(i,j,k)
1137 enddo ; enddo ; enddo
1138 !$OMP parallel do default(shared)
1139 do j=js-1,je ; do k=1,nz ; do i=is,ie
1140 vbt_cor(i,j) = vbt_cor(i,j) + wt_v(i,j,k) * v_cor(i,j,k)
1141 enddo ; enddo ; enddo
1142
1143 ! The gtot arrays are the effective layer-weighted reduced gravities for
1144 ! accelerations across the various faces, with names for the relative
1145 ! locations of the faces to the pressure point. They will have their halos
1146 ! updated later on.
1147 !$OMP parallel do default(shared)
1148 do j=js,je
1149 do k=1,nz ; do i=is-1,ie
1150 gtot_e(i,j) = gtot_e(i,j) + pbce(i,j,k) * wt_u(i,j,k)
1151 gtot_w(i+1,j) = gtot_w(i+1,j) + pbce(i+1,j,k) * wt_u(i,j,k)
1152 enddo ; enddo
1153 enddo
1154 !$OMP parallel do default(shared)
1155 do j=js-1,je
1156 do k=1,nz ; do i=is,ie
1157 gtot_n(i,j) = gtot_n(i,j) + pbce(i,j,k) * wt_v(i,j,k)
1158 gtot_s(i,j+1) = gtot_s(i,j+1) + pbce(i,j+1,k) * wt_v(i,j,k)
1159 enddo ; enddo
1160 enddo
1161
1162 if (cs%BT_OBC%u_OBCs_on_PE) then
1163 do j=js,je ; do i=is-1,ie
1164 if (cs%BT_OBC%u_OBC_type(i,j) > 0) & ! Eastern boundary condition
1165 gtot_w(i+1,j) = gtot_w(i,j) ! Perhaps this should be gtot_E(i,j)?
1166 if (cs%BT_OBC%u_OBC_type(i,j) < 0) & ! Western boundary condition
1167 gtot_e(i,j) = gtot_e(i+1,j) ! Perhaps this should be gtot_W(i+1,j)?
1168 enddo ; enddo
1169 endif
1170 if (cs%BT_OBC%v_OBCs_on_PE) then
1171 do j=js-1,je ; do i=is,ie
1172 if (cs%BT_OBC%v_OBC_type(i,j) > 0) & ! Northern boundary condition
1173 gtot_s(i,j+1) = gtot_s(i,j) !### Should this be gtot_N(i,j) to use wt_v at the same point?
1174 if (cs%BT_OBC%v_OBC_type(i,j) < 0) & ! Southern boundary condition
1175 gtot_n(i,j) = gtot_n(i,j+1) ! Perhaps this should be gtot_S(i,j+1)?
1176 enddo ; enddo
1177 endif
1178
1179 if (cs%calculate_SAL) then
1180 call scalar_sal_sensitivity(cs%SAL_CSp, det_de)
1181 if (cs%tidal_sal_bug) then
1182 dgeo_de = 1.0 + det_de + cs%G_extra
1183 else
1184 dgeo_de = (1.0 - det_de) + cs%G_extra
1185 endif
1186 else
1187 dgeo_de = 1.0 + cs%G_extra
1188 endif
1189
1190 if (nonblock_setup .and. .not.cs%linearized_BT_PV) then
1191 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1192 call complete_group_pass(cs%pass_q_DCor, cs%BT_Domain, clock=id_clock_pass_pre)
1193 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1194 endif
1195
1196 ! Calculate the open areas at the velocity points.
1197 ! The halo updates are needed before Datu is first used, either in set_up_BT_OBC or ubt_Cor.
1198 if (integral_bt_cont) then
1199 call set_local_bt_cont_types(bt_cont, btcl_u, btcl_v, g, us, ms, cs%BT_Domain, 1+ievf-ie, dt_baroclinic=dt)
1200 elseif (use_bt_cont) then
1201 call set_local_bt_cont_types(bt_cont, btcl_u, btcl_v, g, us, ms, cs%BT_Domain, 1+ievf-ie)
1202 else
1203 if (cs%Nonlinear_continuity) then
1204 call find_face_areas(datu, datv, g, gv, us, cs, ms, 1, eta)
1205 else
1206 call find_face_areas(datu, datv, g, gv, us, cs, ms, 1)
1207 endif
1208 endif
1209
1210 ! Set up fields related to the open boundary conditions. These calls include halo updates that
1211 ! must occur on all PEs when there are open boundary conditions anywhere.
1212 if (apply_obcs) then
1213 if (nonblock_setup .and. apply_obc_flather .and. .not.gv%Boussinesq) &
1214 call complete_group_pass(cs%pass_SpV_avg, cs%BT_domain)
1215
1216 dgeo_de_obc = 1.0 ; if (cs%tidal_SAL_Flather) dgeo_de_obc = dgeo_de
1217 call set_up_bt_obc(obc, eta, spv_col_avg, cs%BT_OBC, cs%BT_Domain, g, gv, us, cs, ms, ievf-ie, &
1218 use_bt_cont, integral_bt_cont, dt, datu, datv, btcl_u, btcl_v, dgeo_de_obc)
1219 endif
1220
1221 ! Determine the difference between the sum of the layer fluxes and the
1222 ! barotropic fluxes found from the same input velocities.
1223 if (add_uh0) then
1224 !$OMP parallel do default(shared)
1225 do j=js,je ; do i=is-1,ie ; uhbt(i,j) = 0.0 ; ubt(i,j) = 0.0 ; enddo ; enddo
1226 !$OMP parallel do default(shared)
1227 do j=js-1,je ; do i=is,ie ; vhbt(i,j) = 0.0 ; vbt(i,j) = 0.0 ; enddo ; enddo
1228 if (cs%visc_rem_u_uh0) then
1229 !$OMP parallel do default(shared)
1230 do j=js,je ; do k=1,nz ; do i=is-1,ie
1231 uhbt(i,j) = uhbt(i,j) + uh0(i,j,k)
1232 ubt(i,j) = ubt(i,j) + wt_u(i,j,k) * u_uh0(i,j,k)
1233 enddo ; enddo ; enddo
1234 !$OMP parallel do default(shared)
1235 do j=js-1,je ; do k=1,nz ; do i=is,ie
1236 vhbt(i,j) = vhbt(i,j) + vh0(i,j,k)
1237 vbt(i,j) = vbt(i,j) + wt_v(i,j,k) * v_vh0(i,j,k)
1238 enddo ; enddo ; enddo
1239 else
1240 !$OMP parallel do default(shared)
1241 do j=js,je ; do k=1,nz ; do i=is-1,ie
1242 uhbt(i,j) = uhbt(i,j) + uh0(i,j,k)
1243 ubt(i,j) = ubt(i,j) + cs%frhatu(i,j,k) * u_uh0(i,j,k)
1244 enddo ; enddo ; enddo
1245 !$OMP parallel do default(shared)
1246 do j=js-1,je ; do k=1,nz ; do i=is,ie
1247 vhbt(i,j) = vhbt(i,j) + vh0(i,j,k)
1248 vbt(i,j) = vbt(i,j) + cs%frhatv(i,j,k) * v_vh0(i,j,k)
1249 enddo ; enddo ; enddo
1250 endif
1251 if ((use_bt_cont .or. integral_bt_cont) .and. cs%adjust_BT_cont) then
1252 ! Use the additional input transports to broaden the fits
1253 ! over which the bt_cont_type applies.
1254
1255 ! Fill in the halo data for ubt, vbt, uhbt, and vhbt.
1256 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1257 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1258 call pass_vector(ubt, vbt, cs%BT_Domain, complete=.false., halo=1+ievf-ie)
1259 call pass_vector(uhbt, vhbt, cs%BT_Domain, complete=.true., halo=1+ievf-ie)
1260 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1261 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1262
1263 if (integral_bt_cont) then
1264 call adjust_local_bt_cont_types(ubt, uhbt, vbt, vhbt, btcl_u, btcl_v, &
1265 g, us, ms, 1+ievf-ie, dt_baroclinic=dt)
1266 else
1267 call adjust_local_bt_cont_types(ubt, uhbt, vbt, vhbt, btcl_u, btcl_v, &
1268 g, us, ms, 1+ievf-ie)
1269 endif
1270 endif
1271 if (integral_bt_cont) then
1272 !$OMP parallel do default(shared)
1273 do j=js,je ; do i=is-1,ie
1274 uhbt0(i,j) = uhbt(i,j) - find_uhbt(dt*ubt(i,j), btcl_u(i,j)) * idt
1275 enddo ; enddo
1276 !$OMP parallel do default(shared)
1277 do j=js-1,je ; do i=is,ie
1278 vhbt0(i,j) = vhbt(i,j) - find_vhbt(dt*vbt(i,j), btcl_v(i,j)) * idt
1279 enddo ; enddo
1280 elseif (use_bt_cont) then
1281 !$OMP parallel do default(shared)
1282 do j=js,je ; do i=is-1,ie
1283 uhbt0(i,j) = uhbt(i,j) - find_uhbt(ubt(i,j), btcl_u(i,j))
1284 enddo ; enddo
1285 !$OMP parallel do default(shared)
1286 do j=js-1,je ; do i=is,ie
1287 vhbt0(i,j) = vhbt(i,j) - find_vhbt(vbt(i,j), btcl_v(i,j))
1288 enddo ; enddo
1289 else
1290 !$OMP parallel do default(shared)
1291 do j=js,je ; do i=is-1,ie
1292 uhbt0(i,j) = uhbt(i,j) - datu(i,j)*ubt(i,j)
1293 enddo ; enddo
1294 !$OMP parallel do default(shared)
1295 do j=js-1,je ; do i=is,ie
1296 vhbt0(i,j) = vhbt(i,j) - datv(i,j)*vbt(i,j)
1297 enddo ; enddo
1298 endif
1299 if (cs%BT_OBC%u_OBCs_on_PE) then ! Zero out the reference transport at OBC points
1300 !$OMP parallel do default(shared)
1301 do j=js,je ; do i=is-1,ie ; if (cs%BT_OBC%u_OBC_type(i,j) /= 0) then
1302 uhbt0(i,j) = 0.0
1303 endif ; enddo ; enddo
1304 endif
1305 if (cs%BT_OBC%v_OBCs_on_PE) then !Zero out the reference transport at OBC points
1306 !$OMP parallel do default(shared)
1307 do j=js-1,je ; do i=is,ie ; if (cs%BT_OBC%v_OBC_type(i,j) /= 0) then
1308 vhbt0(i,j) = 0.0
1309 endif ; enddo ; enddo
1310 endif
1311 endif
1312
1313! Calculate the initial barotropic velocities from the layer's velocities.
1314 call btstep_ubt_from_layer(u_in, v_in, wt_u, wt_v, ubt, vbt, g, gv, cs)
1315
1316 uhbt(:,:) = 0.0 ; vhbt(:,:) = 0.0
1317 u_accel_bt(:,:) = 0.0 ; v_accel_bt(:,:) = 0.0
1318
1319 if (apply_obcs .or. (cs%id_ubtdt > 0)) then
1320 do j=js,je ; do i=is-1,ie ; ubt_st(i,j) = ubt(i,j) ; enddo ; enddo
1321 endif
1322 if (apply_obcs .or. (cs%id_vbtdt > 0)) then
1323 do j=js-1,je ; do i=is,ie ; vbt_st(i,j) = vbt(i,j) ; enddo ; enddo
1324 endif
1325
1326! Here the vertical average accelerations due to the Coriolis, advective,
1327! pressure gradient and horizontal viscous terms in the layer momentum
1328! equations are calculated. These will be used to determine the difference
1329! between the accelerations due to the average of the layer equations and the
1330! barotropic calculation.
1331
1332 !$OMP parallel do default(shared)
1333 do j=js,je ; do i=is-1,ie ; if (g%OBCmaskCu(i,j) > 0.0) then
1334 if (cs%nonlin_stress) then
1335 if (gv%Boussinesq) then
1336 htot_avg = 0.5*(max(cs%bathyT(i,j)*gv%Z_to_H + eta(i,j), 0.0) + &
1337 max(cs%bathyT(i+1,j)*gv%Z_to_H + eta(i+1,j), 0.0))
1338 else
1339 htot_avg = 0.5*(eta(i,j) + eta(i+1,j))
1340 endif
1341 if (htot_avg*cs%dy_Cu(i,j) <= 0.0) then
1342 cs%IDatu(i,j) = 0.0
1343 elseif (integral_bt_cont) then
1344 cs%IDatu(i,j) = cs%dy_Cu(i,j) / (max(find_duhbt_du(ubt(i,j)*dt, btcl_u(i,j)), &
1345 cs%dy_Cu(i,j)*htot_avg) )
1346 elseif (use_bt_cont) then ! Reconsider the max and whether there should be some scaling.
1347 cs%IDatu(i,j) = cs%dy_Cu(i,j) / (max(find_duhbt_du(ubt(i,j), btcl_u(i,j)), &
1348 cs%dy_Cu(i,j)*htot_avg) )
1349 else
1350 cs%IDatu(i,j) = 1.0 / htot_avg
1351 endif
1352 endif
1353
1354 bt_force_u(i,j) = forces%taux(i,j) * gv%RZ_to_H * cs%IDatu(i,j)*visc_rem_u(i,j,1)
1355 else
1356 bt_force_u(i,j) = 0.0
1357 endif ; enddo ; enddo
1358 !$OMP parallel do default(shared)
1359 do j=js-1,je ; do i=is,ie ; if (g%OBCmaskCv(i,j) > 0.0) then
1360 if (cs%nonlin_stress) then
1361 if (gv%Boussinesq) then
1362 htot_avg = 0.5*(max(cs%bathyT(i,j)*gv%Z_to_H + eta(i,j), 0.0) + &
1363 max(cs%bathyT(i,j+1)*gv%Z_to_H + eta(i,j+1), 0.0))
1364 else
1365 htot_avg = 0.5*(eta(i,j) + eta(i,j+1))
1366 endif
1367 if (htot_avg*cs%dx_Cv(i,j) <= 0.0) then
1368 cs%IDatv(i,j) = 0.0
1369 elseif (integral_bt_cont) then
1370 cs%IDatv(i,j) = cs%dx_Cv(i,j) / (max(find_dvhbt_dv(vbt(i,j)*dt, btcl_v(i,j)), &
1371 cs%dx_Cv(i,j)*htot_avg) )
1372 elseif (use_bt_cont) then ! Reconsider the max and whether there should be some scaling.
1373 cs%IDatv(i,j) = cs%dx_Cv(i,j) / (max(find_dvhbt_dv(vbt(i,j), btcl_v(i,j)), &
1374 cs%dx_Cv(i,j)*htot_avg) )
1375 else
1376 cs%IDatv(i,j) = 1.0 / htot_avg
1377 endif
1378 endif
1379
1380 bt_force_v(i,j) = forces%tauy(i,j) * gv%RZ_to_H * cs%IDatv(i,j)*visc_rem_v(i,j,1)
1381 else
1382 bt_force_v(i,j) = 0.0
1383 endif ; enddo ; enddo
1384 if (associated(taux_bot) .and. associated(tauy_bot)) then
1385 !$OMP parallel do default(shared)
1386 do j=js,je ; do i=is-1,ie ; if (g%mask2dCu(i,j) > 0.0) then
1387 bt_force_u(i,j) = bt_force_u(i,j) - taux_bot(i,j) * gv%RZ_to_H * cs%IDatu(i,j)
1388 endif ; enddo ; enddo
1389 !$OMP parallel do default(shared)
1390 do j=js-1,je ; do i=is,ie ; if (g%mask2dCv(i,j) > 0.0) then
1391 bt_force_v(i,j) = bt_force_v(i,j) - tauy_bot(i,j) * gv%RZ_to_H * cs%IDatv(i,j)
1392 endif ; enddo ; enddo
1393 endif
1394
1395 ! bc_accel_u & bc_accel_v are only available on the potentially
1396 ! non-symmetric computational domain.
1397 !$OMP parallel do default(shared)
1398 do j=js,je ; do k=1,nz ; do i=isq,ieq
1399 bt_force_u(i,j) = bt_force_u(i,j) + wt_u(i,j,k) * bc_accel_u(i,j,k)
1400 enddo ; enddo ; enddo
1401 !$OMP parallel do default(shared)
1402 do j=jsq,jeq ; do k=1,nz ; do i=is,ie
1403 bt_force_v(i,j) = bt_force_v(i,j) + wt_v(i,j,k) * bc_accel_v(i,j,k)
1404 enddo ; enddo ; enddo
1405
1406 if (cs%gradual_BT_ICs) then
1407 !$OMP parallel do default(shared)
1408 do j=js,je ; do i=is-1,ie
1409 bt_force_u(i,j) = bt_force_u(i,j) + (ubt(i,j) - cs%ubt_IC(i,j)) * idt
1410 ubt(i,j) = cs%ubt_IC(i,j)
1411 if (abs(ubt(i,j)) < cs%vel_underflow) ubt(i,j) = 0.0
1412 enddo ; enddo
1413 !$OMP parallel do default(shared)
1414 do j=js-1,je ; do i=is,ie
1415 bt_force_v(i,j) = bt_force_v(i,j) + (vbt(i,j) - cs%vbt_IC(i,j)) * idt
1416 vbt(i,j) = cs%vbt_IC(i,j)
1417 if (abs(vbt(i,j)) < cs%vel_underflow) vbt(i,j) = 0.0
1418 enddo ; enddo
1419 endif
1420
1421 ! Compute instantaneous tidal velocities and apply frequency-dependent drag.
1422 ! Note that the filtered velocities are only updated during the current predictor step,
1423 ! and are calculated using the barotropic velocity from the previous correction step.
1424 if (cs%use_filter) then
1425 call filt_accum(ubt(g%IsdB:g%IedB,g%jsd:g%jed), ufilt, cs%Time, us, cs%Filt_CS_u)
1426 call filt_accum(vbt(g%isd:g%ied,g%JsdB:g%JedB), vfilt, cs%Time, us, cs%Filt_CS_v)
1427 endif
1428
1429 if (cs%use_filter .and. cs%linear_freq_drag) then
1430 call wave_drag_calc(ufilt, vfilt, drag_u, drag_v, g, cs%Drag_CS)
1431 !$OMP do
1432 do j=js,je ; do i=is-1,ie
1433 htot = 0.5 * (eta(i,j) + eta(i+1,j))
1434 if (gv%Boussinesq) &
1435 htot = htot + 0.5*gv%Z_to_H * (cs%bathyT(i,j) + cs%bathyT(i+1,j))
1436 if (htot > 0.0) then
1437 drag_u(i,j) = drag_u(i,j) / htot
1438 bt_force_u(i,j) = bt_force_u(i,j) - drag_u(i,j)
1439 else
1440 drag_u(i,j) = 0.0
1441 endif
1442 enddo ; enddo
1443 !$OMP do
1444 do j=js-1,je ; do i=is,ie
1445 htot = 0.5 * (eta(i,j) + eta(i,j+1))
1446 if (gv%Boussinesq) &
1447 htot = htot + 0.5*gv%Z_to_H * (cs%bathyT(i,j) + cs%bathyT(i,j+1))
1448 if (htot > 0.0) then
1449 drag_v(i,j) = drag_v(i,j) / htot
1450 bt_force_v(i,j) = bt_force_v(i,j) - drag_v(i,j)
1451 else
1452 drag_v(i,j) = 0.0
1453 endif
1454 enddo ; enddo
1455 endif
1456
1457 ! Mask out the forcing at OBC points
1458 if (cs%BT_OBC%u_OBCs_on_PE) then
1459 !$OMP do
1460 do j=js,je ; do i=is-1,ie
1461 bt_force_u(i,j) = cs%OBCmask_u(i,j) * bt_force_u(i,j)
1462 enddo ; enddo
1463 endif
1464 if (cs%BT_OBC%v_OBCs_on_PE) then
1465 !$OMP do
1466 do j=js-1,je ; do i=is,ie
1467 bt_force_v(i,j) = cs%OBCmask_v(i,j) * bt_force_v(i,j)
1468 enddo ; enddo
1469 endif
1470
1471 if ((isq > is-1) .or. (jsq > js-1)) then
1472 ! Non-symmetric memory is being used, so the edge values need to be
1473 ! filled in with a halo update of a non-symmetric array.
1474 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1475 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1476 tmp_u(:,:) = 0.0 ; tmp_v(:,:) = 0.0
1477 do j=js,je ; do i=isq,ieq ; tmp_u(i,j) = bt_force_u(i,j) ; enddo ; enddo
1478 do j=jsq,jeq ; do i=is,ie ; tmp_v(i,j) = bt_force_v(i,j) ; enddo ; enddo
1479 if (nonblock_setup) then
1480 call start_group_pass(cs%pass_tmp_uv, g%Domain)
1481 else
1482 call do_group_pass(cs%pass_tmp_uv, g%Domain)
1483 do j=jsd,jed ; do i=isdb,iedb ; bt_force_u(i,j) = tmp_u(i,j) ; enddo ; enddo
1484 do j=jsdb,jedb ; do i=isd,ied ; bt_force_v(i,j) = tmp_v(i,j) ; enddo ; enddo
1485 endif
1486 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1487 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1488 endif
1489
1490 if (nonblock_setup) then
1491 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1492 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1493 call start_group_pass(cs%pass_gtot, cs%BT_Domain)
1494 call start_group_pass(cs%pass_ubt_Cor, g%Domain)
1495 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1496 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1497 endif
1498
1499 ! Determine the weighted Coriolis parameters for the neighboring velocities.
1500 call btstep_find_cor(q, dcor_u, dcor_v, f_4_u, f_4_v, isvf, ievf, jsvf, jevf, cs)
1501
1502! Complete the previously initiated message passing.
1503 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1504 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1505 if (nonblock_setup) then
1506 if ((isq > is-1) .or. (jsq > js-1)) then
1507 call complete_group_pass(cs%pass_tmp_uv, g%Domain)
1508 do j=jsd,jed ; do i=isdb,iedb ; bt_force_u(i,j) = tmp_u(i,j) ; enddo ; enddo
1509 do j=jsdb,jedb ; do i=isd,ied ; bt_force_v(i,j) = tmp_v(i,j) ; enddo ; enddo
1510 endif
1511 call complete_group_pass(cs%pass_gtot, cs%BT_Domain)
1512 call complete_group_pass(cs%pass_ubt_Cor, g%Domain)
1513 else
1514 call do_group_pass(cs%pass_gtot, cs%BT_Domain)
1515 call do_group_pass(cs%pass_ubt_Cor, g%Domain)
1516 endif
1517 ! The various elements of gtot are positive definite but directional, so use
1518 ! the polarity arrays to sort out when the directions have shifted.
1519 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf+1
1520 if (cs%ua_polarity(i,j) < 0.0) call swap(gtot_e(i,j), gtot_w(i,j))
1521 if (cs%va_polarity(i,j) < 0.0) call swap(gtot_n(i,j), gtot_s(i,j))
1522 enddo ; enddo
1523
1524 !$OMP parallel do default(shared)
1525 do j=js,je ; do i=is-1,ie
1526 cor_ref_u(i,j) = &
1527 (((f_4_u(4,i,j) * vbt_cor(i+1,j)) + (f_4_u(1,i,j) * vbt_cor(i ,j-1))) + &
1528 ((f_4_u(3,i,j) * vbt_cor(i ,j)) + (f_4_u(2,i,j) * vbt_cor(i+1,j-1))))
1529 enddo ; enddo
1530 !$OMP parallel do default(shared)
1531 do j=js-1,je ; do i=is,ie
1532 cor_ref_v(i,j) = -1.0 * &
1533 (((f_4_v(1,i,j) * ubt_cor(i-1,j)) + (f_4_v(4,i,j) * ubt_cor(i ,j+1))) + &
1534 ((f_4_v(2,i,j) * ubt_cor(i ,j)) + (f_4_v(3,i,j) * ubt_cor(i-1,j+1))))
1535 enddo ; enddo
1536
1537 ! Now start new halo updates.
1538 if (nonblock_setup) then
1539 if (.not.use_bt_cont) &
1540 call start_group_pass(cs%pass_Dat_uv, cs%BT_Domain)
1541
1542 ! The following halo update is not needed without wide halos. RWH
1543 call start_group_pass(cs%pass_force_hbt0_Cor_ref, cs%BT_Domain)
1544 endif
1545 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1546 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1547 !$OMP parallel default(shared) private(u_max_cor,uint_cor,v_max_cor,vint_cor,eta_cor_max,Htot)
1548 !$OMP do
1549 do j=js-1,je+1 ; do i=is-1,ie ; av_rem_u(i,j) = 0.0 ; enddo ; enddo
1550 !$OMP do
1551 do j=js-1,je ; do i=is-1,ie+1 ; av_rem_v(i,j) = 0.0 ; enddo ; enddo
1552 !$OMP do
1553 do j=js,je ; do k=1,nz ; do i=is-1,ie
1554 av_rem_u(i,j) = av_rem_u(i,j) + cs%frhatu(i,j,k) * visc_rem_u(i,j,k)
1555 enddo ; enddo ; enddo
1556 !$OMP do
1557 do j=js-1,je ; do k=1,nz ; do i=is,ie
1558 av_rem_v(i,j) = av_rem_v(i,j) + cs%frhatv(i,j,k) * visc_rem_v(i,j,k)
1559 enddo ; enddo ; enddo
1560 if (cs%strong_drag) then
1561 !$OMP do
1562 do j=js,je ; do i=is-1,ie
1563 bt_rem_u(i,j) = g%mask2dCu(i,j) * &
1564 ((nstep * av_rem_u(i,j)) / (1.0 + (nstep-1)*av_rem_u(i,j)))
1565 enddo ; enddo
1566 !$OMP do
1567 do j=js-1,je ; do i=is,ie
1568 bt_rem_v(i,j) = g%mask2dCv(i,j) * &
1569 ((nstep * av_rem_v(i,j)) / (1.0 + (nstep-1)*av_rem_v(i,j)))
1570 enddo ; enddo
1571 else
1572 !$OMP do
1573 do j=js,je ; do i=is-1,ie
1574 bt_rem_u(i,j) = 0.0
1575 if (g%mask2dCu(i,j) * av_rem_u(i,j) > 0.0) &
1576 bt_rem_u(i,j) = g%mask2dCu(i,j) * (av_rem_u(i,j)**instep)
1577 enddo ; enddo
1578 !$OMP do
1579 do j=js-1,je ; do i=is,ie
1580 bt_rem_v(i,j) = 0.0
1581 if (g%mask2dCv(i,j) * av_rem_v(i,j) > 0.0) &
1582 bt_rem_v(i,j) = g%mask2dCv(i,j) * (av_rem_v(i,j)**instep)
1583 enddo ; enddo
1584 endif
1585 if (cs%linear_wave_drag) then
1586 !$OMP do
1587 do j=js,je ; do i=is-1,ie ; if (g%mask2dCu(i,j) * cs%lin_drag_u(i,j) > 0.0) then
1588 htot = 0.5 * (eta(i,j) + eta(i+1,j))
1589 if (gv%Boussinesq) &
1590 htot = htot + 0.5*gv%Z_to_H * (cs%bathyT(i,j) + cs%bathyT(i+1,j))
1591 ! If Htot==0., linear wave drag is not used and Rayleigh_u = 0.0 (from initialization)
1592 ! and bt_rem_u is unmodified.
1593 if (htot > 0.0) then
1594 bt_rem_u(i,j) = bt_rem_u(i,j) * (htot / (htot + cs%lin_drag_u(i,j) * dtbt))
1595 rayleigh_u(i,j) = cs%lin_drag_u(i,j) / htot
1596 endif
1597 endif ; enddo ; enddo
1598 !$OMP do
1599 do j=js-1,je ; do i=is,ie ; if (g%mask2dCv(i,j) * cs%lin_drag_v(i,j) > 0.0) then
1600 htot = 0.5 * (eta(i,j) + eta(i,j+1))
1601 if (gv%Boussinesq) &
1602 htot = htot + 0.5*gv%Z_to_H * (cs%bathyT(i,j) + cs%bathyT(i,j+1))
1603 ! If Htot==0., linear wave drag is not used and Rayleigh_v = 0.0 (from initialization)
1604 ! and bt_rem_v is unmodified.
1605 if (htot > 0.0) then
1606 bt_rem_v(i,j) = bt_rem_v(i,j) * (htot / (htot + cs%lin_drag_v(i,j) * dtbt))
1607 rayleigh_v(i,j) = cs%lin_drag_v(i,j) / htot
1608 endif
1609 endif ; enddo ; enddo
1610 endif
1611
1612 ! Avoid changing the velocities at OBC points due to non-OBC calculations.
1613 if (cs%BT_OBC%u_OBCs_on_PE) then
1614 !$OMP do
1615 do j=js,je ; do i=is-1,ie ; if (cs%BT_OBC%u_OBC_type(i,j) /= 0) then
1616 bt_rem_u(i,j) = 1.0
1617 endif ; enddo ; enddo
1618 endif
1619 if (cs%BT_OBC%v_OBCs_on_PE) then
1620 !$OMP do
1621 do j=js-1,je ; do i=is,ie ; if (cs%BT_OBC%v_OBC_type(i,j) /= 0) then
1622 bt_rem_v(i,j) = 1.0
1623 endif ; enddo ; enddo
1624 endif
1625
1626 ! Set the mass source, after first initializing the halos to 0.
1627 !$OMP do
1628 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf+1 ; eta_src(i,j) = 0.0 ; enddo ; enddo
1629 if (cs%bound_BT_corr) then ; if ((use_bt_cont.or.integral_bt_cont) .and. cs%BT_cont_bounds) then
1630 do j=js,je ; do i=is,ie ; if (g%mask2dT(i,j) > 0.0) then
1631 if (cs%eta_cor(i,j) > 0.0) then
1632 ! Limit the source (outward) correction to be a fraction the mass that
1633 ! can be transported out of the cell by velocities with a CFL number of CFL_cor.
1634 if (integral_bt_cont) then
1635 uint_cor = g%dxT(i,j) * cs%maxCFL_BT_cont
1636 vint_cor = g%dyT(i,j) * cs%maxCFL_BT_cont
1637 eta_cor_max = (cs%IareaT(i,j) * &
1638 (((find_uhbt(uint_cor, btcl_u(i,j)) + dt*uhbt0(i,j)) - &
1639 (find_uhbt(-uint_cor, btcl_u(i-1,j)) + dt*uhbt0(i-1,j))) + &
1640 ((find_vhbt(vint_cor, btcl_v(i,j)) + dt*vhbt0(i,j)) - &
1641 (find_vhbt(-vint_cor, btcl_v(i,j-1)) + dt*vhbt0(i,j-1))) ))
1642 else ! (use_BT_Cont) then
1643 u_max_cor = g%dxT(i,j) * (cs%maxCFL_BT_cont*idt)
1644 v_max_cor = g%dyT(i,j) * (cs%maxCFL_BT_cont*idt)
1645 eta_cor_max = dt * (cs%IareaT(i,j) * &
1646 (((find_uhbt(u_max_cor, btcl_u(i,j)) + uhbt0(i,j)) - &
1647 (find_uhbt(-u_max_cor, btcl_u(i-1,j)) + uhbt0(i-1,j))) + &
1648 ((find_vhbt(v_max_cor, btcl_v(i,j)) + vhbt0(i,j)) - &
1649 (find_vhbt(-v_max_cor, btcl_v(i,j-1)) + vhbt0(i,j-1))) ))
1650 endif
1651 cs%eta_cor(i,j) = min(cs%eta_cor(i,j), max(0.0, eta_cor_max))
1652 else
1653 ! Limit the sink (inward) correction to the amount of mass that is already inside the cell.
1654 htot = eta(i,j)
1655 if (gv%Boussinesq) htot = cs%bathyT(i,j)*gv%Z_to_H + eta(i,j)
1656
1657 cs%eta_cor(i,j) = max(cs%eta_cor(i,j), -max(0.0,htot))
1658 endif
1659 endif ; enddo ; enddo
1660 else ; do j=js,je ; do i=is,ie
1661 if (abs(cs%eta_cor(i,j)) > dt*cs%eta_cor_bound(i,j)) &
1662 cs%eta_cor(i,j) = sign(dt*cs%eta_cor_bound(i,j), cs%eta_cor(i,j))
1663 enddo ; enddo ; endif ; endif
1664 !$OMP do
1665 do j=js,je ; do i=is,ie
1666 eta_src(i,j) = g%mask2dT(i,j) * (instep * cs%eta_cor(i,j))
1667 enddo ; enddo
1668 !$OMP end parallel
1669
1670 if (cs%dynamic_psurf) then
1671 ice_is_rigid = (associated(forces%rigidity_ice_u) .and. &
1672 associated(forces%rigidity_ice_v))
1673 h_min_dyn = cs%Dmin_dyn_psurf
1674 if (ice_is_rigid .and. use_bt_cont) &
1675 call bt_cont_to_face_areas(bt_cont, datu, datv, g, us, ms, halo=0)
1676 if (ice_is_rigid) then
1677 if (gv%Boussinesq) then
1678 h_to_z = gv%H_to_Z
1679 else
1680 h_to_z = gv%H_to_RZ / cs%Rho_BT_lin
1681 endif
1682 !$OMP parallel do default(shared) private(Idt_max2,H_eff_dx2,dyn_coef_max,ice_strength)
1683 do j=js,je ; do i=is,ie
1684 ! First determine the maximum stable value for dyn_coef_eta.
1685
1686 ! This estimate of the maximum stable time step is pretty accurate for
1687 ! gravity waves, but it is a conservative estimate since it ignores the
1688 ! stabilizing effect of the bottom drag.
1689 idt_max2 = 0.5 * (dgeo_de * (1.0 + 2.0*cs%bebt)) * (g%IareaT(i,j) * &
1690 (((gtot_e(i,j) * (datu(i,j)*g%IdxCu(i,j))) + &
1691 (gtot_w(i,j) * (datu(i-1,j)*g%IdxCu(i-1,j)))) + &
1692 ((gtot_n(i,j) * (datv(i,j)*g%IdyCv(i,j))) + &
1693 (gtot_s(i,j) * (datv(i,j-1)*g%IdyCv(i,j-1))))) + &
1694 ((g%Coriolis2Bu(i,j) + g%Coriolis2Bu(i-1,j-1)) + &
1695 (g%Coriolis2Bu(i-1,j) + g%Coriolis2Bu(i,j-1))) * cs%BT_Coriolis_scale**2 )
1696 h_eff_dx2 = max(h_min_dyn * ((g%IdxT(i,j)**2) + (g%IdyT(i,j)**2)), &
1697 g%IareaT(i,j) * &
1698 (((datu(i,j)*g%IdxCu(i,j)) + (datu(i-1,j)*g%IdxCu(i-1,j))) + &
1699 ((datv(i,j)*g%IdyCv(i,j)) + (datv(i,j-1)*g%IdyCv(i,j-1))) ) )
1700 dyn_coef_max = cs%const_dyn_psurf * max(0.0, 1.0 - dtbt**2 * idt_max2) / &
1701 (dtbt**2 * h_eff_dx2)
1702
1703 ! ice_strength has units of [L2 Z-1 T-2 ~> m s-2]. rigidity_ice_[uv] has units of [L4 Z-1 T-1 ~> m3 s-1].
1704 ice_strength = ((forces%rigidity_ice_u(i,j) + forces%rigidity_ice_u(i-1,j)) + &
1705 (forces%rigidity_ice_v(i,j) + forces%rigidity_ice_v(i,j-1))) / &
1706 (cs%ice_strength_length**2 * dtbt)
1707
1708 ! Units of dyn_coef: [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1]
1709 dyn_coef_eta(i,j) = min(dyn_coef_max, ice_strength * h_to_z)
1710 enddo ; enddo ; endif
1711 endif
1712
1713 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1714 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1715 if (nonblock_setup) then
1716 call start_group_pass(cs%pass_eta_bt_rem, cs%BT_Domain)
1717 ! The following halo update is not needed without wide halos. RWH
1718 else
1719 call do_group_pass(cs%pass_eta_bt_rem, cs%BT_Domain)
1720 if (.not.use_bt_cont) &
1721 call do_group_pass(cs%pass_Dat_uv, cs%BT_Domain)
1722 call do_group_pass(cs%pass_force_hbt0_Cor_ref, cs%BT_Domain)
1723 endif
1724 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1725 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1726
1727 ! Complete all of the outstanding halo updates.
1728 if (nonblock_setup) then
1729 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1730 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
1731
1732 if (.not.use_bt_cont) call complete_group_pass(cs%pass_Dat_uv, cs%BT_Domain)
1733 call complete_group_pass(cs%pass_force_hbt0_Cor_ref, cs%BT_Domain)
1734 call complete_group_pass(cs%pass_eta_bt_rem, cs%BT_Domain)
1735
1736 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
1737 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
1738 endif
1739
1740 if (cs%debug) then
1741 call uvchksum("BT [uv]hbt", uhbt, vhbt, cs%debug_BT_HI, haloshift=0, &
1742 unscale=us%s_to_T*us%L_to_m**2*gv%H_to_m)
1743 call uvchksum("BT Initial [uv]bt", ubt, vbt, cs%debug_BT_HI, haloshift=0, unscale=us%L_T_to_m_s)
1744 call hchksum(eta, "BT Initial eta", cs%debug_BT_HI, haloshift=0, unscale=gv%H_to_MKS)
1745 call uvchksum("BT BT_force_[uv]", bt_force_u, bt_force_v, &
1746 cs%debug_BT_HI, haloshift=0, unscale=us%L_T2_to_m_s2)
1747 if (interp_eta_pf) then
1748 call hchksum(eta_pf_1, "BT eta_PF_1", cs%debug_BT_HI, haloshift=0, unscale=gv%H_to_MKS)
1749 call hchksum(d_eta_pf, "BT d_eta_PF", cs%debug_BT_HI, haloshift=0, unscale=gv%H_to_MKS)
1750 else
1751 call hchksum(eta_pf, "BT eta_PF", cs%debug_BT_HI, haloshift=0, unscale=gv%H_to_MKS)
1752 call hchksum(eta_pf_in, "BT eta_PF_in", g%HI,haloshift=0, unscale=gv%H_to_MKS)
1753 endif
1754 if (cs%linearized_BT_PV) then
1755 call bchksum(cs%q_D, "BT PV (q_D)", cs%debug_BT_HI, haloshift=0, symmetric=.true., unscale=us%s_to_T/gv%H_to_MKS)
1756 else
1757 call bchksum(q, "BT PV (q)", cs%debug_BT_HI, haloshift=0, symmetric=.true., unscale=us%s_to_T/gv%H_to_MKS)
1758 endif
1759 call uvchksum("BT DCor_[uv]", dcor_u, dcor_v, g%HI, haloshift=0, &
1760 symmetric=.true., omit_corners=.true., scalar_pair=.true., unscale=gv%H_to_MKS)
1761 call uvchksum("BT Cor_ref_[uv]", cor_ref_u, cor_ref_v, cs%debug_BT_HI, haloshift=0, unscale=us%L_T2_to_m_s2)
1762 call uvchksum("BT [uv]hbt0", uhbt0, vhbt0, cs%debug_BT_HI, haloshift=0, &
1763 unscale=us%L_to_m**2*us%s_to_T*gv%H_to_m)
1764 if (.not. use_bt_cont) then
1765 call uvchksum("BT Dat[uv]", datu, datv, cs%debug_BT_HI, haloshift=1, unscale=us%L_to_m*gv%H_to_m)
1766 endif
1767 call uvchksum("BT wt_[uv]", wt_u, wt_v, g%HI, haloshift=0, &
1768 symmetric=.true., omit_corners=.true., scalar_pair=.true.)
1769 call uvchksum("BT frhat[uv]", cs%frhatu, cs%frhatv, g%HI, haloshift=0, &
1770 symmetric=.true., omit_corners=.true., scalar_pair=.true.)
1771 call uvchksum("BT visc_rem_[uv]", visc_rem_u, visc_rem_v, g%HI, haloshift=0, &
1772 symmetric=.true., omit_corners=.true., scalar_pair=.true.)
1773 call uvchksum("BT bc_accel_[uv]", bc_accel_u, bc_accel_v, g%HI, haloshift=0, unscale=us%L_T2_to_m_s2)
1774 call uvchksum("BT IDat[uv]", cs%IDatu, cs%IDatv, g%HI, haloshift=0, &
1775 unscale=gv%m_to_H, scalar_pair=.true.)
1776 call uvchksum("BT visc_rem_[uv]", visc_rem_u, visc_rem_v, g%HI, &
1777 haloshift=1, scalar_pair=.true.)
1778
1779 if (apply_obcs) then
1780 call uvchksum("BT_OBC%[uv]bt_outer", cs%BT_OBC%ubt_outer, cs%BT_OBC%vbt_outer, cs%debug_BT_HI, &
1781 symmetric=.true., omit_corners=.true., unscale=us%L_T_to_m_s)
1782 if (allocated(cs%BT_OBC%SSH_outer_u) .and. allocated(cs%BT_OBC%SSH_outer_v)) &
1783 call uvchksum("BT_OBC%SSH_outer[uv]", cs%BT_OBC%SSH_outer_u, cs%BT_OBC%SSH_outer_v, cs%debug_BT_HI, &
1784 symmetric=.true., omit_corners=.true., unscale=us%Z_to_m, scalar_pair=.true.)
1785 if (allocated(cs%BT_OBC%Cg_u) .and. allocated(cs%BT_OBC%Cg_v)) &
1786 call uvchksum("BT_OBC%Cg_[uv]", cs%BT_OBC%Cg_u, cs%BT_OBC%Cg_v, cs%debug_BT_HI, &
1787 symmetric=.true., omit_corners=.true., unscale=us%L_T_to_m_s, scalar_pair=.true.)
1788 if (allocated(cs%BT_OBC%dZ_u) .and. allocated(cs%BT_OBC%dZ_v)) &
1789 call uvchksum("BT_OBC%dZ_[uv]", cs%BT_OBC%dZ_u, cs%BT_OBC%dZ_v, cs%debug_BT_HI, &
1790 symmetric=.true., omit_corners=.true., unscale=us%Z_to_m, scalar_pair=.true.)
1791 endif
1792 endif
1793
1794 if (query_averaging_enabled(cs%diag)) then
1795 if (cs%id_eta_st > 0) call post_data(cs%id_eta_st, eta(isd:ied,jsd:jed), cs%diag)
1796 if (cs%id_ubt_st > 0) call post_data(cs%id_ubt_st, ubt(isdb:iedb,jsd:jed), cs%diag)
1797 if (cs%id_vbt_st > 0) call post_data(cs%id_vbt_st, vbt(isd:ied,jsdb:jedb), cs%diag)
1798 endif
1799
1800 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
1801 if (id_clock_calc > 0) call cpu_clock_begin(id_clock_calc)
1802
1803 if (cs%dt_bt_filter >= 0.0) then
1804 dt_filt = 0.5 * max(0.0, min(cs%dt_bt_filter, 2.0*dt))
1805 else
1806 dt_filt = 0.5 * max(0.0, dt * min(-cs%dt_bt_filter, 2.0))
1807 endif
1808 nfilter = ceiling(dt_filt / dtbt)
1809
1810 if ( nstep+nfilter<=0 ) call mom_error(fatal, &
1811 "btstep: number of barotropic step (nstep+nfilter) is 0")
1812 if ( cs%bt_limit_integral_transport .and. nstep-nfilter<=0 ) call mom_error(fatal, &
1813 "btstep: barotropic filter steps too large (nstep-nfilter) is 0")
1814
1815 ! Set up the normalized weights for the filtered velocity.
1816 sum_wt_vel = 0.0 ; sum_wt_eta = 0.0 ; sum_wt_accel = 0.0 ; sum_wt_trans = 0.0
1817 allocate(wt_vel(nstep+nfilter)) ; allocate(wt_eta(nstep+nfilter))
1818 allocate(wt_trans(nstep+nfilter+1)) ; allocate(wt_accel(nstep+nfilter+1))
1819 allocate(wt_accel2(nstep+nfilter+1))
1820 do n=1,nstep+nfilter
1821 ! Modify this to use a different filter...
1822
1823 ! This is a filter that ramps down linearly over a time dt_filt.
1824 if ( (n==nstep) .or. (dt_filt - abs(n-nstep)*dtbt >= 0.0)) then
1825 wt_vel(n) = 1.0 ; wt_eta(n) = 1.0
1826 elseif (dtbt + dt_filt - abs(n-nstep)*dtbt > 0.0) then
1827 wt_vel(n) = 1.0 + (dt_filt / dtbt) - abs(n-nstep) ; wt_eta(n) = wt_vel(n)
1828 else
1829 wt_vel(n) = 0.0 ; wt_eta(n) = 0.0
1830 endif
1831 ! This is a simple stepfunction filter.
1832 ! if (n < nstep-nfilter) then ; wt_vel(n) = 0.0 ; else ; wt_vel(n) = 1.0 ; endif
1833 ! wt_eta(n) = wt_vel(n)
1834
1835 ! The rest should not be changed.
1836 sum_wt_vel = sum_wt_vel + wt_vel(n) ; sum_wt_eta = sum_wt_eta + wt_eta(n)
1837 enddo
1838 wt_trans(nstep+nfilter+1) = 0.0 ; wt_accel(nstep+nfilter+1) = 0.0
1839 do n=nstep+nfilter,1,-1
1840 wt_trans(n) = wt_trans(n+1) + wt_eta(n)
1841 wt_accel(n) = wt_accel(n+1) + wt_vel(n)
1842 sum_wt_accel = sum_wt_accel + wt_accel(n) ; sum_wt_trans = sum_wt_trans + wt_trans(n)
1843 enddo
1844 ! Normalize the weights.
1845 i_sum_wt_vel = 1.0 / sum_wt_vel ; i_sum_wt_accel = 1.0 / sum_wt_accel
1846 i_sum_wt_eta = 1.0 / sum_wt_eta ; i_sum_wt_trans = 1.0 / sum_wt_trans
1847 do n=1,nstep+nfilter
1848 wt_vel(n) = wt_vel(n) * i_sum_wt_vel
1849 if (cs%answer_date < 20190101) then
1850 wt_accel2(n) = wt_accel(n)
1851 ! wt_trans(n) = wt_trans(n) * I_sum_wt_trans
1852 else
1853 wt_accel2(n) = wt_accel(n) * i_sum_wt_accel
1854 wt_trans(n) = wt_trans(n) * i_sum_wt_trans
1855 endif
1856 wt_accel(n) = wt_accel(n) * i_sum_wt_accel
1857 wt_eta(n) = wt_eta(n) * i_sum_wt_eta
1858 enddo
1859 if (cs%answer_date < 20190101) then
1860 ! Recalculate the sum of the weights even that they may have been renormalized already.
1861 sum_wt_vel = 0.0 ; sum_wt_eta = 0.0 ; sum_wt_trans = 0.0 ; sum_wt_accel = 0.0
1862 do n=1,nstep+nfilter
1863 sum_wt_vel = sum_wt_vel + wt_vel(n)
1864 sum_wt_eta = sum_wt_eta + wt_eta(n)
1865 sum_wt_accel = sum_wt_accel + wt_accel2(n)
1866 sum_wt_trans = sum_wt_trans + wt_trans(n)
1867 enddo
1868 i_sum_wt_vel = 1.0 / sum_wt_vel ; i_sum_wt_eta = 1.0 / sum_wt_eta
1869 i_sum_wt_accel = 1.0 / sum_wt_accel ; i_sum_wt_trans = 1.0 / sum_wt_trans
1870 else
1871 i_sum_wt_vel = 1.0 ; i_sum_wt_eta = 1.0 ; i_sum_wt_accel = 1.0 ; i_sum_wt_trans = 1.0
1872 endif
1873
1874 ! March the barotropic solver through all of its time steps.
1875 call btstep_timeloop(eta, ubt, vbt, uhbt0, datu, btcl_u, vhbt0, datv, btcl_v, eta_ic, &
1876 eta_pf_1, d_eta_pf, eta_src, dyn_coef_eta, uhbtav, vhbtav, u_accel_bt, v_accel_bt, &
1877 f_4_u, f_4_v, bt_rem_u, bt_rem_v, &
1878 bt_force_u, bt_force_v, cor_ref_u, cor_ref_v, rayleigh_u, rayleigh_v, &
1879 eta_pf, gtot_e, gtot_w, gtot_n, gtot_s, spv_col_avg, dgeo_de, &
1880 eta_sum, eta_wtd, ubt_wtd, vbt_wtd, coru_avg, pfu_avg, ldu_avg, corv_avg, pfv_avg, &
1881 ldv_avg, use_bt_cont, interp_eta_pf, find_etaav, dt, dtbt, nstep, nfilter, &
1882 wt_vel, wt_eta, wt_accel, wt_trans, wt_accel2, adp, cs%BT_OBC, cs, g, ms, gv, us)
1883
1884 if (id_clock_calc > 0) call cpu_clock_end(id_clock_calc)
1885 if (id_clock_calc_post > 0) call cpu_clock_begin(id_clock_calc_post)
1886
1887 if (find_etaav) then ; do j=js,je ; do i=is,ie
1888 etaav(i,j) = eta_sum(i,j) * i_sum_wt_accel
1889 enddo ; enddo ; endif
1890 do j=js-1,je+1 ; do i=is-1,ie+1 ; e_anom(i,j) = 0.0 ; enddo ; enddo
1891 if (interp_eta_pf) then
1892 do j=js,je ; do i=is,ie
1893 e_anom(i,j) = dgeo_de * (0.5 * (eta(i,j) + eta_in(i,j)) - &
1894 (eta_pf_1(i,j) + 0.5*d_eta_pf(i,j)))
1895 enddo ; enddo
1896 else
1897 do j=js,je ; do i=is,ie
1898 e_anom(i,j) = dgeo_de * (0.5 * (eta(i,j) + eta_in(i,j)) - eta_pf(i,j))
1899 enddo ; enddo
1900 endif
1901 if (apply_obcs) then
1902 ! This block of code may be unnecessary because e_anom is only used for accelerations that
1903 ! are then recalculated at OBC points.
1904 if (cs%BT_OBC%u_OBCs_on_PE) then ! copy back the value for u-points on the boundary.
1905 !GOMP parallel do default(shared)
1906 do j=js,je ; do i=is-1,ie
1907 if (cs%BT_OBC%u_OBC_type(i,j) > 0) e_anom(i+1,j) = e_anom(i,j) ! OBC_DIRECTION_E
1908 if (cs%BT_OBC%u_OBC_type(i,j) < 0) e_anom(i,j) = e_anom(i+1,j) ! OBC_DIRECTION_W
1909 enddo ; enddo
1910 endif
1911
1912 if (cs%BT_OBC%v_OBCs_on_PE) then ! copy back the value for v-points on the boundary.
1913 !GOMP parallel do default(shared)
1914 do j=js-1,je ; do i=is,ie
1915 if (cs%BT_OBC%v_OBC_type(i,j) > 0) e_anom(i,j+1) = e_anom(i,j) ! OBC_DIRECTION_N
1916 if (cs%BT_OBC%v_OBC_type(i,j) < 0) e_anom(i,j) = e_anom(i,j+1) ! OBC_DIRECTION_S
1917 enddo ; enddo
1918 endif
1919 endif
1920
1921 ! Note that it is possible that eta_out and eta_in are the same array.
1922 do j=js,je ; do i=is,ie
1923 eta_out(i,j) = eta_wtd(i,j) * i_sum_wt_eta
1924 enddo ; enddo
1925
1926 ! Accumulator is updated at the end of every baroclinic time step.
1927 ! Harmonic analysis will not be performed of a field that is not registered.
1928 if (associated(cs%HA_CSp) .and. find_etaav) then
1929 call ha_accum('ubt', ubt, cs%Time, g, cs%HA_CSp)
1930 call ha_accum('vbt', vbt, cs%Time, g, cs%HA_CSp)
1931 endif
1932
1933 if (id_clock_calc_post > 0) call cpu_clock_end(id_clock_calc_post)
1934 if (id_clock_pass_post > 0) call cpu_clock_begin(id_clock_pass_post)
1935 if (g%nonblocking_updates) then
1936 call start_group_pass(cs%pass_e_anom, g%Domain)
1937 else
1938 if (find_etaav) call do_group_pass(cs%pass_etaav, g%Domain)
1939 call do_group_pass(cs%pass_e_anom, g%Domain)
1940 endif
1941 if (id_clock_pass_post > 0) call cpu_clock_end(id_clock_pass_post)
1942 if (id_clock_calc_post > 0) call cpu_clock_begin(id_clock_calc_post)
1943
1944 ! Find or store the weighted time-mean velocities and transports.
1945 if (cs%answer_date < 20190101) then
1946 do j=js,je ; do i=is-1,ie
1947 cs%ubtav(i,j) = cs%ubtav(i,j) * i_sum_wt_trans
1948 uhbtav(i,j) = uhbtav(i,j) * i_sum_wt_trans
1949 ubt_wtd(i,j) = ubt_wtd(i,j) * i_sum_wt_vel
1950 enddo ; enddo
1951
1952 do j=js-1,je ; do i=is,ie
1953 cs%vbtav(i,j) = cs%vbtav(i,j) * i_sum_wt_trans
1954 vhbtav(i,j) = vhbtav(i,j) * i_sum_wt_trans
1955 vbt_wtd(i,j) = vbt_wtd(i,j) * i_sum_wt_vel
1956 enddo ; enddo
1957 endif
1958
1959 if (cs%use_filter .and. cs%linear_freq_drag) then ! Apply frequency-dependent drag
1960 !$OMP do
1961 do j=js,je ; do i=is-1,ie
1962 u_accel_bt(i,j) = u_accel_bt(i,j) - drag_u(i,j)
1963 enddo ; enddo
1964 !$OMP do
1965 do j=js-1,je ; do i=is,ie
1966 v_accel_bt(i,j) = v_accel_bt(i,j) - drag_v(i,j)
1967 enddo ; enddo
1968
1969 if ((cs%id_LDu_bt > 0) .or. (associated(adp%bt_lwd_u))) then ; do j=js,je ; do i=is-1,ie
1970 ldu_avg(i,j) = ldu_avg(i,j) - drag_u(i,j)
1971 enddo ; enddo ; endif
1972 if ((cs%id_LDv_bt > 0) .or. (associated(adp%bt_lwd_v))) then ; do j=js-1,je ; do i=is,ie
1973 ldv_avg(i,j) = ldv_avg(i,j) - drag_v(i,j)
1974 enddo ; enddo ; endif
1975 endif
1976
1977 if (id_clock_calc_post > 0) call cpu_clock_end(id_clock_calc_post)
1978 if (id_clock_pass_post > 0) call cpu_clock_begin(id_clock_pass_post)
1979 if (g%nonblocking_updates) then
1980 call complete_group_pass(cs%pass_e_anom, g%Domain)
1981 if (find_etaav) call start_group_pass(cs%pass_etaav, g%Domain)
1982 call start_group_pass(cs%pass_ubta_uhbta, g%DoMain)
1983 else
1984 call do_group_pass(cs%pass_ubta_uhbta, g%Domain)
1985 endif
1986 if (id_clock_pass_post > 0) call cpu_clock_end(id_clock_pass_post)
1987 if (id_clock_calc_post > 0) call cpu_clock_begin(id_clock_calc_post)
1988
1989 if (cs%strong_drag .and. cs%rescale_strong_drag) then
1990 do j=js,je ; do i=is-1,ie
1991 if (g%mask2dCu(i,j) * av_rem_u(i,j) > 0.0) &
1992 u_accel_bt(i,j) = u_accel_bt(i,j) * min(bt_rem_u(i,j)**nstep / av_rem_u(i,j), 1.0)
1993 enddo ; enddo
1994 do j=js-1,je ; do i=is,ie
1995 if (g%mask2dCv(i,j) * av_rem_v(i,j) > 0.0) &
1996 v_accel_bt(i,j) = v_accel_bt(i,j) * min(bt_rem_v(i,j)**nstep / av_rem_v(i,j), 1.0)
1997 enddo ; enddo
1998 endif
1999
2000 ! Now calculate each layer's accelerations.
2001 call btstep_layer_accel(dt, u_accel_bt, v_accel_bt, pbce, gtot_e, gtot_w, gtot_n, gtot_s, &
2002 e_anom, g, gv, cs, accel_layer_u, accel_layer_v)
2003
2004 if (apply_obcs) then
2005 ! Correct the accelerations at OBC velocity points, but only in the
2006 ! symmetric-memory computational domain, not in the wide halo regions.
2007 if (cs%BT_OBC%u_OBCs_on_PE) then ; do j=js,je ; do i=is-1,ie
2008 if (cs%BT_OBC%u_OBC_type(i,j) /= 0) then
2009 u_accel_bt(i,j) = (ubt_wtd(i,j) - ubt_st(i,j)) / dt
2010 do k=1,nz ; accel_layer_u(i,j,k) = u_accel_bt(i,j) ; enddo
2011 endif
2012 enddo ; enddo ; endif
2013 if (cs%BT_OBC%v_OBCs_on_PE) then ; do j=js-1,je ; do i=is,ie
2014 if (cs%BT_OBC%v_OBC_type(i,j) /= 0) then
2015 v_accel_bt(i,j) = (vbt_wtd(i,j) - vbt_st(i,j)) / dt
2016 do k=1,nz ; accel_layer_v(i,j,k) = v_accel_bt(i,j) ; enddo
2017 endif
2018 enddo ; enddo ; endif
2019 endif
2020
2021 if (id_clock_calc_post > 0) call cpu_clock_end(id_clock_calc_post)
2022
2023 ! Calculate diagnostic quantities.
2024 if (query_averaging_enabled(cs%diag)) then
2025
2026 if (cs%gradual_BT_ICs) then
2027 do j=js,je ; do i=is-1,ie ; cs%ubt_IC(i,j) = ubt_wtd(i,j) ; enddo ; enddo
2028 do j=js-1,je ; do i=is,ie ; cs%vbt_IC(i,j) = vbt_wtd(i,j) ; enddo ; enddo
2029 endif
2030
2031 ! Calculate various time-averaged barotropic diagnostics.
2032 if (cs%answer_date >= 20190101) then
2033 if (cs%id_PFu_bt > 0) call post_data(cs%id_PFu_bt, pfu_avg, cs%diag)
2034 if (cs%id_PFv_bt > 0) call post_data(cs%id_PFv_bt, pfv_avg, cs%diag)
2035 if (cs%id_Coru_bt > 0) call post_data(cs%id_Coru_bt, coru_avg, cs%diag)
2036 if (cs%id_Corv_bt > 0) call post_data(cs%id_Corv_bt, corv_avg, cs%diag)
2037 if (cs%id_LDu_bt > 0) call post_data(cs%id_LDu_bt, ldu_avg, cs%diag)
2038 if (cs%id_LDv_bt > 0) call post_data(cs%id_LDv_bt, ldv_avg, cs%diag)
2039 else ! if (CS%answer_date < 20190101) then
2040 if (cs%id_PFu_bt > 0) then
2041 do j=js,je ; do i=is-1,ie
2042 pfu_avg(i,j) = pfu_avg(i,j) * i_sum_wt_accel
2043 enddo ; enddo
2044 call post_data(cs%id_PFu_bt, pfu_avg, cs%diag)
2045 endif
2046 if (cs%id_PFv_bt > 0) then
2047 do j=js-1,je ; do i=is,ie
2048 pfv_avg(i,j) = pfv_avg(i,j) * i_sum_wt_accel
2049 enddo ; enddo
2050 call post_data(cs%id_PFv_bt, pfv_avg, cs%diag)
2051 endif
2052 if (cs%id_Coru_bt > 0) then
2053 do j=js,je ; do i=is-1,ie
2054 coru_avg(i,j) = coru_avg(i,j) * i_sum_wt_accel
2055 enddo ; enddo
2056 call post_data(cs%id_Coru_bt, coru_avg, cs%diag)
2057 endif
2058 if (cs%id_Corv_bt > 0) then
2059 do j=js-1,je ; do i=is,ie
2060 corv_avg(i,j) = corv_avg(i,j) * i_sum_wt_accel
2061 enddo ; enddo
2062 call post_data(cs%id_Corv_bt, corv_avg, cs%diag)
2063 endif
2064 endif
2065
2066 ! Diagnostics for time tendency
2067 if (cs%id_ubtdt > 0) then
2068 do j=js,je ; do i=is-1,ie
2069 ubt_dt(i,j) = (ubt_wtd(i,j) - ubt_st(i,j))*idt
2070 enddo ; enddo
2071 call post_data(cs%id_ubtdt, ubt_dt(isdb:iedb,jsd:jed), cs%diag)
2072 endif
2073 if (cs%id_vbtdt > 0) then
2074 do j=js-1,je ; do i=is,ie
2075 vbt_dt(i,j) = (vbt_wtd(i,j) - vbt_st(i,j))*idt
2076 enddo ; enddo
2077 call post_data(cs%id_vbtdt, vbt_dt(isd:ied,jsdb:jedb), cs%diag)
2078 endif
2079
2080 ! Copy decomposed barotropic accelerations to ADp
2081 if (associated(adp%bt_pgf_u)) then
2082 ! Note that CS%IdxCu is 0 at OBC points, so ADp%bt_pgf_u is zeroed out there.
2083 do k=1,nz ; do j=js,je ; do i=is-1,ie
2084 adp%bt_pgf_u(i,j,k) = pfu_avg(i,j) - &
2085 (((pbce(i+1,j,k) - gtot_w(i+1,j)) * e_anom(i+1,j)) - &
2086 ((pbce(i,j,k) - gtot_e(i,j)) * e_anom(i,j))) * cs%IdxCu(i,j)
2087 enddo ; enddo ; enddo
2088 endif
2089 if (associated(adp%bt_pgf_v)) then
2090 ! Note that CS%IdyCv is 0 at OBC points, so ADp%bt_pgf_v is zeroed out there.
2091 do k=1,nz ; do j=js-1,je ; do i=is,ie
2092 adp%bt_pgf_v(i,j,k) = pfv_avg(i,j) - &
2093 (((pbce(i,j+1,k) - gtot_s(i,j+1)) * e_anom(i,j+1)) - &
2094 ((pbce(i,j,k) - gtot_n(i,j)) * e_anom(i,j))) * cs%IdyCv(i,j)
2095 enddo ; enddo ; enddo
2096 endif
2097
2098 if (associated(adp%bt_cor_u)) then ; do j=js,je ; do i=is-1,ie
2099 adp%bt_cor_u(i,j) = coru_avg(i,j)
2100 enddo ; enddo ; endif
2101 if (associated(adp%bt_cor_v)) then ; do j=js-1,je ; do i=is,ie
2102 adp%bt_cor_v(i,j) = corv_avg(i,j)
2103 enddo ; enddo ; endif
2104
2105 if (associated(adp%bt_lwd_u)) then ; do j=js,je ; do i=is-1,ie
2106 adp%bt_lwd_u(i,j) = ldu_avg(i,j)
2107 enddo ; enddo ; endif
2108 if (associated(adp%bt_lwd_v)) then ; do j=js-1,je ; do i=is,ie
2109 adp%bt_lwd_v(i,j) = ldv_avg(i,j)
2110 enddo ; enddo ; endif
2111
2112 if (cs%id_ubtforce > 0) call post_data(cs%id_ubtforce, bt_force_u(isdb:iedb,jsd:jed), cs%diag)
2113 if (cs%id_vbtforce > 0) call post_data(cs%id_vbtforce, bt_force_v(isd:ied,jsdb:jedb), cs%diag)
2114 if (cs%id_uaccel > 0) call post_data(cs%id_uaccel, u_accel_bt(isdb:iedb,jsd:jed), cs%diag)
2115 if (cs%id_vaccel > 0) call post_data(cs%id_vaccel, v_accel_bt(isd:ied,jsdb:jedb), cs%diag)
2116
2117 if (cs%id_eta_cor > 0) call post_data(cs%id_eta_cor, cs%eta_cor, cs%diag)
2118 if (cs%id_eta_bt > 0) call post_data(cs%id_eta_bt, eta_out, cs%diag) ! - G%Z_ref?
2119 if (cs%id_gtotn > 0) call post_data(cs%id_gtotn, gtot_n(isd:ied,jsd:jed), cs%diag)
2120 if (cs%id_gtots > 0) call post_data(cs%id_gtots, gtot_s(isd:ied,jsd:jed), cs%diag)
2121 if (cs%id_gtote > 0) call post_data(cs%id_gtote, gtot_e(isd:ied,jsd:jed), cs%diag)
2122 if (cs%id_gtotw > 0) call post_data(cs%id_gtotw, gtot_w(isd:ied,jsd:jed), cs%diag)
2123 if (cs%id_ubt > 0) call post_data(cs%id_ubt, ubt_wtd(isdb:iedb,jsd:jed), cs%diag)
2124 if (cs%id_vbt > 0) call post_data(cs%id_vbt, vbt_wtd(isd:ied,jsdb:jedb), cs%diag)
2125 if (cs%id_ubtav > 0) call post_data(cs%id_ubtav, cs%ubtav, cs%diag)
2126 if (cs%id_vbtav > 0) call post_data(cs%id_vbtav, cs%vbtav, cs%diag)
2127 if (cs%id_visc_rem_u > 0) call post_data(cs%id_visc_rem_u, visc_rem_u, cs%diag)
2128 if (cs%id_visc_rem_v > 0) call post_data(cs%id_visc_rem_v, visc_rem_v, cs%diag)
2129
2130 if (cs%id_frhatu > 0) call post_data(cs%id_frhatu, cs%frhatu, cs%diag)
2131 if (cs%id_uhbt > 0) call post_data(cs%id_uhbt, uhbtav, cs%diag)
2132 if (cs%id_frhatv > 0) call post_data(cs%id_frhatv, cs%frhatv, cs%diag)
2133 if (cs%id_vhbt > 0) call post_data(cs%id_vhbt, vhbtav, cs%diag)
2134 if (cs%id_uhbt0 > 0) call post_data(cs%id_uhbt0, uhbt0(isdb:iedb,jsd:jed), cs%diag)
2135 if (cs%id_vhbt0 > 0) call post_data(cs%id_vhbt0, vhbt0(isd:ied,jsdb:jedb), cs%diag)
2136
2137 if (cs%id_frhatu1 > 0) call post_data(cs%id_frhatu1, cs%frhatu1, cs%diag)
2138 if (cs%id_frhatv1 > 0) call post_data(cs%id_frhatv1, cs%frhatv1, cs%diag)
2139
2140 if (cs%id_bt_rem_u > 0) call post_data(cs%id_bt_rem_u, bt_rem_u(isdb:iedb,jsd:jed), cs%diag)
2141 if (cs%id_bt_rem_v > 0) call post_data(cs%id_bt_rem_v, bt_rem_v(isd:ied,jsdb:jedb), cs%diag)
2142 if (cs%id_etaPF_anom > 0) call post_data(cs%id_etaPF_anom, e_anom(isd:ied,jsd:jed), cs%diag)
2143
2144 if (use_bt_cont) then
2145 if (cs%id_BTC_FA_u_EE > 0) call post_data(cs%id_BTC_FA_u_EE, bt_cont%FA_u_EE, cs%diag)
2146 if (cs%id_BTC_FA_u_E0 > 0) call post_data(cs%id_BTC_FA_u_E0, bt_cont%FA_u_E0, cs%diag)
2147 if (cs%id_BTC_FA_u_W0 > 0) call post_data(cs%id_BTC_FA_u_W0, bt_cont%FA_u_W0, cs%diag)
2148 if (cs%id_BTC_FA_u_WW > 0) call post_data(cs%id_BTC_FA_u_WW, bt_cont%FA_u_WW, cs%diag)
2149 if (cs%id_BTC_uBT_EE > 0) call post_data(cs%id_BTC_uBT_EE, bt_cont%uBT_EE, cs%diag)
2150 if (cs%id_BTC_uBT_WW > 0) call post_data(cs%id_BTC_uBT_WW, bt_cont%uBT_WW, cs%diag)
2151 if (cs%id_BTC_FA_u_rat0 > 0) then
2152 tmp_u(:,:) = 0.0
2153 do j=js,je ; do i=is-1,ie
2154 if ((g%mask2dCu(i,j) > 0.0) .and. (bt_cont%FA_u_W0(i,j) > 0.0)) then
2155 tmp_u(i,j) = (bt_cont%FA_u_E0(i,j)/ bt_cont%FA_u_W0(i,j))
2156 else
2157 tmp_u(i,j) = 1.0
2158 endif
2159 enddo ; enddo
2160 call post_data(cs%id_BTC_FA_u_rat0, tmp_u, cs%diag)
2161 endif
2162 if (cs%id_BTC_FA_v_NN > 0) call post_data(cs%id_BTC_FA_v_NN, bt_cont%FA_v_NN, cs%diag)
2163 if (cs%id_BTC_FA_v_N0 > 0) call post_data(cs%id_BTC_FA_v_N0, bt_cont%FA_v_N0, cs%diag)
2164 if (cs%id_BTC_FA_v_S0 > 0) call post_data(cs%id_BTC_FA_v_S0, bt_cont%FA_v_S0, cs%diag)
2165 if (cs%id_BTC_FA_v_SS > 0) call post_data(cs%id_BTC_FA_v_SS, bt_cont%FA_v_SS, cs%diag)
2166 if (cs%id_BTC_vBT_NN > 0) call post_data(cs%id_BTC_vBT_NN, bt_cont%vBT_NN, cs%diag)
2167 if (cs%id_BTC_vBT_SS > 0) call post_data(cs%id_BTC_vBT_SS, bt_cont%vBT_SS, cs%diag)
2168 if (cs%id_BTC_FA_v_rat0 > 0) then
2169 tmp_v(:,:) = 0.0
2170 do j=js-1,je ; do i=is,ie
2171 if ((g%mask2dCv(i,j) > 0.0) .and. (bt_cont%FA_v_S0(i,j) > 0.0)) then
2172 tmp_v(i,j) = (bt_cont%FA_v_N0(i,j)/ bt_cont%FA_v_S0(i,j))
2173 else
2174 tmp_v(i,j) = 1.0
2175 endif
2176 enddo ; enddo
2177 call post_data(cs%id_BTC_FA_v_rat0, tmp_v, cs%diag)
2178 endif
2179 if (cs%id_BTC_FA_h_rat0 > 0) then
2180 tmp_h(:,:) = 0.0
2181 do j=js,je ; do i=is,ie
2182 tmp_h(i,j) = 1.0
2183 if ((g%mask2dCu(i,j) > 0.0) .and. (bt_cont%FA_u_W0(i,j) > 0.0) .and. (bt_cont%FA_u_E0(i,j) > 0.0)) then
2184 if (bt_cont%FA_u_W0(i,j) > bt_cont%FA_u_E0(i,j)) then
2185 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_u_W0(i,j)/ bt_cont%FA_u_E0(i,j)))
2186 else
2187 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_u_E0(i,j)/ bt_cont%FA_u_W0(i,j)))
2188 endif
2189 endif
2190 if ((g%mask2dCu(i-1,j) > 0.0) .and. (bt_cont%FA_u_W0(i-1,j) > 0.0) .and. (bt_cont%FA_u_E0(i-1,j) > 0.0)) then
2191 if (bt_cont%FA_u_W0(i-1,j) > bt_cont%FA_u_E0(i-1,j)) then
2192 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_u_W0(i-1,j)/ bt_cont%FA_u_E0(i-1,j)))
2193 else
2194 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_u_E0(i-1,j)/ bt_cont%FA_u_W0(i-1,j)))
2195 endif
2196 endif
2197 if ((g%mask2dCv(i,j) > 0.0) .and. (bt_cont%FA_v_S0(i,j) > 0.0) .and. (bt_cont%FA_v_N0(i,j) > 0.0)) then
2198 if (bt_cont%FA_v_S0(i,j) > bt_cont%FA_v_N0(i,j)) then
2199 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_v_S0(i,j)/ bt_cont%FA_v_N0(i,j)))
2200 else
2201 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_v_N0(i,j)/ bt_cont%FA_v_S0(i,j)))
2202 endif
2203 endif
2204 if ((g%mask2dCv(i,j-1) > 0.0) .and. (bt_cont%FA_v_S0(i,j-1) > 0.0) .and. (bt_cont%FA_v_N0(i,j-1) > 0.0)) then
2205 if (bt_cont%FA_v_S0(i,j-1) > bt_cont%FA_v_N0(i,j-1)) then
2206 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_v_S0(i,j-1)/ bt_cont%FA_v_N0(i,j-1)))
2207 else
2208 tmp_h(i,j) = max(tmp_h(i,j), (bt_cont%FA_v_N0(i,j-1)/ bt_cont%FA_v_S0(i,j-1)))
2209 endif
2210 endif
2211 enddo ; enddo
2212 call post_data(cs%id_BTC_FA_h_rat0, tmp_h, cs%diag)
2213 endif
2214 endif
2215
2216 if (cs%id_SSH_u_OBC > 0) call post_data(cs%id_SSH_u_OBC, cs%BT_OBC%SSH_outer_u(isdb:iedb,jsd:jed), cs%diag)
2217 if (cs%id_SSH_v_OBC > 0) call post_data(cs%id_SSH_v_OBC, cs%BT_OBC%SSH_outer_v(isd:ied,jsdb:jedb), cs%diag)
2218 if (cs%id_ubt_OBC > 0) call post_data(cs%id_ubt_OBC, cs%BT_OBC%ubt_outer(isdb:iedb,jsd:jed), cs%diag)
2219 if (cs%id_vbt_OBC > 0) call post_data(cs%id_vbt_OBC, cs%BT_OBC%vbt_outer(isd:ied,jsdb:jedb), cs%diag)
2220 else
2221 if (cs%id_frhatu1 > 0) cs%frhatu1(:,:,:) = cs%frhatu(:,:,:)
2222 if (cs%id_frhatv1 > 0) cs%frhatv1(:,:,:) = cs%frhatv(:,:,:)
2223 endif
2224
2225 if (associated(adp%diag_hfrac_u)) then
2226 do k=1,nz ; do j=js,je ; do i=is-1,ie
2227 adp%diag_hfrac_u(i,j,k) = cs%frhatu(i,j,k)
2228 enddo ; enddo ; enddo
2229 endif
2230 if (associated(adp%diag_hfrac_v)) then
2231 do k=1,nz ; do j=js-1,je ; do i=is,ie
2232 adp%diag_hfrac_v(i,j,k) = cs%frhatv(i,j,k)
2233 enddo ; enddo ; enddo
2234 endif
2235
2236 if (use_bt_cont .and. associated(adp%diag_hu)) then
2237 do k=1,nz ; do j=js,je ; do i=is-1,ie
2238 adp%diag_hu(i,j,k) = bt_cont%h_u(i,j,k)
2239 enddo ; enddo ; enddo
2240 endif
2241 if (use_bt_cont .and. associated(adp%diag_hv)) then
2242 do k=1,nz ; do j=js-1,je ; do i=is,ie
2243 adp%diag_hv(i,j,k) = bt_cont%h_v(i,j,k)
2244 enddo ; enddo ; enddo
2245 endif
2246 if (associated(adp%visc_rem_u)) then
2247 do k=1,nz ; do j=js,je ; do i=is-1,ie
2248 adp%visc_rem_u(i,j,k) = visc_rem_u(i,j,k)
2249 enddo ; enddo ; enddo
2250 endif
2251 if (associated(adp%visc_rem_v)) then
2252 do k=1,nz ; do j=js-1,je ; do i=is,ie
2253 adp%visc_rem_v(i,j,k) = visc_rem_v(i,j,k)
2254 enddo ; enddo ; enddo
2255 endif
2256
2257 if (g%nonblocking_updates) then
2258 if (find_etaav) call complete_group_pass(cs%pass_etaav, g%Domain)
2259 call complete_group_pass(cs%pass_ubta_uhbta, g%Domain)
2260 endif
2261
2262 deallocate(wt_vel, wt_eta, wt_trans, wt_accel, wt_accel2)
2263
2264end subroutine btstep
2265
2266!> Update the barotropic solver through multiple time steps.
2267subroutine btstep_timeloop(eta, ubt, vbt, uhbt0, Datu, BTCL_u, vhbt0, Datv, BTCL_v, eta_IC, &
2268 eta_PF_1, d_eta_PF, eta_src, dyn_coef_eta, uhbtav, vhbtav, u_accel_bt, v_accel_bt, &
2269 f_4_u, f_4_v, bt_rem_u, bt_rem_v, &
2270 BT_force_u, BT_force_v, Cor_ref_u, Cor_ref_v, Rayleigh_u, Rayleigh_v, &
2271 eta_PF, gtot_E, gtot_W, gtot_N, gtot_S, SpV_col_avg, dgeo_de, &
2272 eta_sum, eta_wtd, ubt_wtd, vbt_wtd, Coru_avg, PFu_avg, LDu_avg, Corv_avg, PFv_avg, &
2273 LDv_avg, use_BT_cont, interp_eta_PF, find_etaav, dt, dtbt, nstep, nfilter, &
2274 wt_vel, wt_eta, wt_accel, wt_trans, wt_accel2, ADp, BT_OBC, CS, G, MS, GV, US)
2275
2276 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
2277 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure (inout to allow for halo updates)
2278 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of
2279 !! the argument arrays.
2280 real, dimension(SZIW_(CS),SZJW_(CS)), target, intent(inout) :: &
2281 eta !< The barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2]
2282 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
2283 ubt !< The zonal barotropic velocity [L T-1 ~> m s-1]
2284 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
2285 vbt !< The meridional barotropic velocity [L T-1 ~> m s-1]
2286 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
2287 uhbt0 !< The difference between the sum of the layer zonal thickness flux and the
2288 !! barotropic thickness flux using the same velocity [H L2 T-1 ~> m3 s-1 or kg s-1]
2289 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
2290 Datu !< Basin depth at u-velocity grid points times the y-grid spacing [H L ~> m2 or kg m-1]
2291 type(local_bt_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: &
2292 BTCL_u !< Structure of information used for a dynamic estimate of the face areas at u-points.
2293 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2294 vhbt0 !< The difference between the sum of the layer meridional thickness flux and the
2295 !! barotropic thickness flux using the same velocity [H L2 T-1 ~> m3 s-1 or kg s-1]
2296 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
2297 Datv !< Basin depth at v-velocity grid points times the x-grid spacing [H L ~> m2 or kg m-1]
2298 type(local_bt_cont_v_type), dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: &
2299 BTCL_v !< Structure of information used for a dynamic estimate of the face areas at v-points
2300 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2301 eta_IC !< A local copy of the initial 2-D eta field (eta_in) [H ~> m or kg m-2]
2302 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2303 eta_PF_1 !< The initial value of eta_PF, when interp_eta_PF is true [H ~> m or kg m-2]
2304 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2305 d_eta_PF !< The change in eta_PF over the barotropic time stepping when
2306 !! interp_eta_PF is true [H ~> m or kg m-2]
2307 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2308 eta_src !< The source of eta per barotropic timestep [H ~> m or kg m-2]
2309 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2310 dyn_coef_eta !< The coefficient relating the changes in eta to the dynamic surface pressure
2311 !! under rigid ice [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
2312 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: &
2313 uhbtav !< the barotropic zonal volume or mass fluxes averaged through the barotropic
2314 !! steps [H L2 T-1 ~> m3 s-1 or kg s-1].
2315 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: &
2316 vhbtav !< the barotropic meridional volume or mass fluxes averaged through the barotropic
2317 !! steps [H L2 T-1 ~> m3 s-1 or kg s-1].
2318 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
2319 u_accel_bt !! The difference between the zonal acceleration from the
2320 !< barotropic calculation and BT_force_v [L T-2 ~> m s-2].
2321 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
2322 v_accel_bt !< The difference between the meridional acceleration from the
2323 !! barotropic calculation and BT_force_v [L T-2 ~> m s-2].
2324 real, dimension(4,SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
2325 f_4_u !< The terms giving the contribution to the Coriolis acceleration at a zonal
2326 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
2327 !! These are the products of thicknesses at v points and appropriately staggered
2328 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
2329 !! they are approximately f / 4. The 4 values on the innermost loop are for
2330 !! v-velocities to the southwest, southeast, northwest and northeast.
2331 real, dimension(4,SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2332 f_4_v !< The terms giving the contribution to the Coriolis acceleration at a meridional
2333 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
2334 !! These are the products of thicknesses at u points and appropriately staggered
2335 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
2336 !! they are approximately f / 4. The 4 values on the innermost loop are for
2337 !! u-velocities to the southwest, southeast, northwest and northeast.
2338 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
2339 bt_rem_u !< The fraction of the barotropic zonal velocity that remains after a time step,
2340 !! the rest being lost to bottom drag [nondim]. bt_rem_v is between 0 and 1.
2341 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2342 bt_rem_v !< The fraction of the barotropic meridional velocity that remains after a time step,
2343 !! the rest being lost to bottom drag [nondim]. bt_rem_v is between 0 and 1.
2344 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
2345 BT_force_u !< The vertical average of all of the v-accelerations that are
2346 !! not explicitly included in the barotropic equation [L T-2 ~> m s-2]
2347 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2348 BT_force_v !< The vertical average of all of the v-accelerations that are
2349 !! not explicitly included in the barotropic equation [L T-2 ~> m s-2]
2350 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
2351 Cor_ref_u !< The meridional barotropic Coriolis acceleration due
2352 !! to the reference velocities [L T-2 ~> m s-2].
2353 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2354 Cor_ref_v !< The meridional barotropic Coriolis acceleration due
2355 !! to the reference velocities [L T-2 ~> m s-2].
2356 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
2357 Rayleigh_u !< A Rayleigh drag timescale operating at u-points for drag parameterizations
2358 !! that introduced directly into the barotropic solver rather than coming
2359 !! in via the visc_rem_u arrays from the layered equations [T-1 ~> s-1]
2360 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
2361 Rayleigh_v !< A Rayleigh drag timescale operating at v-points for drag parameterizations
2362 !! that introduced directly into the barotropic solver rather than coming
2363 !! in via the visc_rem_v arrays from the layered equations [T-1 ~> s-1]
2364 real, dimension(SZIW_(CS),SZJW_(CS)), intent(inout) :: &
2365 eta_PF !< The 2-D eta field (either SSH anomaly or column mass anomaly) that was used to
2366 !! calculate the input pressure gradient accelerations [H ~> m or kg m-2]
2367 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2368 gtot_E !< The effective total reduced gravity used to relate free surface height
2369 !! deviations to pressure forces (including GFS and baroclinic contributions)
2370 !! in the barotropic momentum equations half a grid-point to the east of a
2371 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
2372 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2373 gtot_W !< The effective total reduced gravity used to relate free surface height
2374 !! deviations to pressure forces (including GFS and baroclinic contributions)
2375 !! in the barotropic momentum equations half a grid-point to the west of a
2376 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2]
2377 !! (See Hallberg, J Comp Phys 1997 for a discussion of gtot_E and gtot_W.)
2378 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2379 gtot_N !< The effective total reduced gravity used to relate free surface height
2380 !! deviations to pressure forces (including GFS and baroclinic contributions)
2381 !! in the barotropic momentum equations half a grid-point to the north of a
2382 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2]
2383 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
2384 gtot_S !< The effective total reduced gravity used to relate free surface height
2385 !! deviations to pressure forces (including GFS and baroclinic contributions)
2386 !! in the barotropic momentum equations half a grid-point to the south of a
2387 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2]
2388 !! (See Hallberg, J Comp Phys 1997 for a discussion of gtot_E and gtot_W.)
2389 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: &
2390 SpV_col_avg !< The column average specific volume [R-1 ~> m3 kg-1]
2391 real, intent(in) :: dgeo_de !< The constant of proportionality between geopotential and
2392 !! sea surface height [nondim]. It is of order 1, but for stability this
2393 !! may be made larger than the physical problem would suggest.
2394 real, dimension(SZIW_(CS),SZJW_(CS)), intent(out) :: &
2395 eta_sum !< eta summed across the timesteps [H ~> m or kg m-2]
2396 real, dimension(SZIW_(CS),SZJW_(CS)), intent(out) :: &
2397 eta_wtd !< A weighted estimate used to calculate eta_out [H ~> m or kg m-2]
2398 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: &
2399 ubt_wtd !< A weighted sum used to find the filtered final ubt [L T-1 ~> m s-1]
2400 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: &
2401 vbt_wtd !< A weighted sum used to find the filtered final vbt [L T-1 ~> m s-1]
2402 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: &
2403 Coru_avg !< The average zonal barotropic Coriolis acceleration [L T-2 ~> m s-2]
2404 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: &
2405 PFu_avg !< The average zonal barotropic pressure gradient force [L T-2 ~> m s-2]
2406 real, dimension(SZIB_(G),SZJ_(G)), intent(out) :: &
2407 LDu_avg !< The average zonal barotropic linear wave drag acceleration [L T-2 ~> m s-2]
2408 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: &
2409 Corv_avg !< The average meridional barotropic Coriolis acceleration [L T-2 ~> m s-2]
2410 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: &
2411 PFv_avg !< The average meridional barotropic pressure gradient force [L T-2 ~> m s-2]
2412 real, dimension(SZI_(G),SZJB_(G)), intent(out) :: &
2413 LDv_avg !< The average meridional barotropic linear wave drag acceleration [L T-2 ~> m s-2]
2414 logical, intent(in) :: use_BT_cont !< If true, use the information in the bt_cont_types to
2415 !! calculate the mass transports
2416 logical, intent(in) :: interp_eta_PF !< If true, interpolate the reference value of eta used
2417 !! to calculate the pressure force with time.
2418 logical, intent(in) :: find_etaav !< If true, diagnose the time mean value of eta
2419 real, intent(in) :: dt !< The time increment to integrate over [T ~> s]
2420 real, intent(in) :: dtbt !< The barotropic time step [T ~> s]
2421 integer, intent(in) :: nstep !< The number of barotropic time steps to take to cover the specified time interval
2422 integer, intent(in) :: nfilter !< The number of extra barotropic steps to take to allow for time filtering
2423 real, dimension(nstep+nfilter), intent(in) :: &
2424 wt_vel !< The raw or relative weights of each of the barotropic timesteps
2425 !! in determining the average velocities [nondim]
2426 real, dimension(nstep+nfilter), intent(in) :: &
2427 wt_eta !< The raw or relative weights of each of the barotropic timesteps
2428 !! in determining the average eta [nondim]
2429 real, dimension(nstep+nfilter+1), intent(in) :: &
2430 wt_accel !< The raw or relative weights of each of the barotropic timesteps
2431 !! in determining the average accelerations [nondim]
2432 real, dimension(nstep+nfilter+1), intent(in) :: &
2433 wt_trans !< The raw or relative weights of each of the barotropic timesteps
2434 !! in determining the average transports [nondim]
2435 real, dimension(nstep+nfilter+1), intent(in) :: &
2436 wt_accel2 !< Potentially un-normalized relative weights of each of the
2437 !! barotropic timesteps in determining the average accelerations [nondim]
2438 type(accel_diag_ptrs), pointer :: ADp !< Acceleration diagnostic pointers
2439 type(bt_obc_type), intent(in) :: BT_OBC !< A structure with the private barotropic arrays
2440 !! related to the open boundary conditions,
2441 !! with time evolving data stored via set_up_BT_OBC
2442 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
2443 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
2444
2445 ! Local variables
2446 real, dimension(SZIBW_(CS),SZJW_(CS)) :: &
2447 uhbt, & ! The zonal barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1]
2448 ubt_prev, & ! The starting value of ubt in a barotropic step [L T-1 ~> m s-1]
2449 ubt_trans, & ! The latest value of ubt used for a transport [L T-1 ~> m s-1]
2450 PFu, & ! The zonal pressure force acceleration [L T-2 ~> m s-2]
2451 Cor_u, & ! The zonal Coriolis acceleration [L T-2 ~> m s-2]
2452 ubt_int, & ! The running time integral of ubt over the time steps [L ~> m]
2453 uhbt_int, & ! The running time integral of uhbt over the time steps [H L2 ~> m3 or kg]
2454 ubt_int_prev, & ! Previous value of time-integrated velocity stored for OBCs [L ~> m]
2455 uhbt_int_prev ! Previous value of time-integrated transport stored for integral_BT_cont [H L2 ~> m3 or kg]
2456 real, dimension(SZIW_(CS),SZJBW_(CS)) :: &
2457 vhbt, & ! The meridional barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1]
2458 vbt_prev, & ! The starting value of vbt in a barotropic step [L T-1 ~> m s-1]
2459 vbt_trans, & ! The latest value of vbt used for a transport [L T-1 ~> m s-1]
2460 PFv, & ! The meridional pressure force acceleration [L T-2 ~> m s-2]
2461 Cor_v, & ! The meridional Coriolis acceleration [L T-2 ~> m s-2]
2462 vbt_int, & ! The running time integral of vbt over the time steps [L ~> m]
2463 vhbt_int, & ! The running time integral of vhbt over the time steps [H L2 ~> m3 or kg]
2464 vbt_int_prev, & ! Previous value of time-integrated velocity stored for OBCs [L ~> m]
2465 vhbt_int_prev ! Previous value of time-integrated transport stored for integral_BT_cont [H L2 ~> m3 or kg]
2466 real, target, dimension(SZIW_(CS),SZJW_(CS)) :: &
2467 eta_pred ! A predictor value of eta [H ~> m or kg m-2] like eta
2468 real, dimension(SZIW_(CS),SZJW_(CS)) :: &
2469 p_surf_dyn, & !< A dynamic surface pressure under rigid ice [L2 T-2 ~> m2 s-2]
2470 cfl_ltd_vol !< The volume available after removing sinks used to limit uhbt_int and vhbt_int [H L2 ~> m3 or kg]
2471 real, dimension(SZI_(G),SZJ_(G)) :: &
2472 eta_anom_PF ! The eta anomalies used to find the pressure force anomalies [H ~> m or kg m-2]
2473 real :: wt_end ! The weighting of the final value of eta_PF [nondim]
2474 real :: Instep ! The inverse of the number of barotropic time steps to take [nondim]
2475 real :: trans_wt1, trans_wt2 ! The weights used to compute ubt_trans and vbt_trans [nondim]
2476 real :: Idtbt ! The inverse of the barotropic time step [T-1 ~> s-1]
2477 type(time_type) :: &
2478 time_bt_start, & ! The starting time of the barotropic steps.
2479 time_step_end, & ! The end time of a barotropic step.
2480 time_end_in ! The end time for diagnostics when this routine started.
2481 real :: dtbt_diag ! The nominal barotropic time step used in hifreq diagnostics [T ~> s]
2482 ! dtbt_diag = dt/(nstep+nfilter)
2483 real :: time_int_in ! The diagnostics' time interval when this routine started [s]
2484 real :: be_proj ! The fractional amount by which velocities are projected
2485 ! when project_velocity is true [nondim]. For now be_proj is set
2486 ! to equal bebt, as they have similar roles and meanings.
2487 real :: eta_cor_multiplier ! Increases the rate of applying CS%eta_cor so that the mass
2488 ! source is all used up by the beginning of the filtering [nondim]
2489 real :: eta_acc ! Change due to divergence of mass transport [H ~> m or kg m-2]
2490 logical :: do_hifreq_output ! If true, output occurs every barotropic step.
2491 logical :: do_ave ! If true, diagnostics are enabled on this step.
2492 logical :: evolving_face_areas
2493 logical :: v_first ! If true, update the v-velocity first with the present loop iteration
2494 logical :: integral_BT_cont ! If true, update the barotropic continuity equation directly
2495 ! from the initial condition using the time-integrated barotropic velocity.
2496 character(len=200) :: mesg
2497 integer :: isv, iev, jsv, jev ! The valid array size at the end of a step.
2498 integer :: isvf, ievf, jsvf, jevf ! The fullest range of array indices that could be used.
2499 integer :: num_cycles ! The number of timesteps before a halo update is needed.
2500 integer :: stencil ! The stencil size of the algorithm, often 1 or 2.
2501 integer :: err_count ! A counter to limit the volume of error messages written to stdout.
2502 integer :: i, j, n, is, ie, js, je
2503 integer :: debug_halo ! The halo size to use for debugging checksums
2504 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
2505
2506 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
2507 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
2508 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
2509 err_count = 0
2510
2511 ! Figure out the fullest arrays that could be updated.
2512 stencil = max(1, cs%min_stencil)
2513 if ((.not.use_bt_cont) .and. cs%Nonlinear_continuity .and. (cs%Nonlin_cont_update_period > 0)) &
2514 stencil = max(2, cs%min_stencil)
2515
2516 num_cycles = 1
2517 if (cs%use_wide_halos) &
2518 num_cycles = min((is-cs%isdw) / stencil, (js-cs%jsdw) / stencil)
2519 isvf = is - (num_cycles-1)*stencil ; ievf = ie + (num_cycles-1)*stencil
2520 jsvf = js - (num_cycles-1)*stencil ; jevf = je + (num_cycles-1)*stencil
2521
2522 integral_bt_cont = use_bt_cont .and. cs%integral_BT_cont
2523 evolving_face_areas = (.not.use_bt_cont) .and. cs%Nonlinear_continuity .and. &
2524 (cs%Nonlin_cont_update_period > 0)
2525 instep = 1.0 / real(nstep)
2526 idtbt = 1.0 / dtbt
2527
2528 !--- setup the weight when computing vbt_trans and ubt_trans
2529 if (cs%BT_project_velocity) then
2530 be_proj = cs%bebt
2531 trans_wt1 = (1.0 + be_proj) ; trans_wt2 = -be_proj
2532 else
2533 trans_wt1 = cs%bebt ; trans_wt2 = (1.0-cs%bebt)
2534 endif
2535
2536 ! Manage diagnostics
2537 do_ave = query_averaging_enabled(cs%diag) .and. &
2538 ((cs%id_PFu_bt > 0) .or. (cs%id_Coru_bt > 0) .or. (cs%id_LDu_bt > 0) .or. &
2539 (cs%id_PFv_bt > 0) .or. (cs%id_Corv_bt > 0) .or. (cs%id_LDv_bt > 0) .or. &
2540 associated(adp%bt_pgf_u) .or. associated(adp%bt_cor_u) .or. associated(adp%bt_lwd_u) .or. &
2541 associated(adp%bt_pgf_v) .or. associated(adp%bt_cor_v) .or. associated(adp%bt_lwd_v))
2542
2543 do_hifreq_output = .false.
2544 if ((cs%id_ubt_hifreq > 0) .or. (cs%id_vbt_hifreq > 0) .or. &
2545 (cs%id_eta_hifreq > 0) .or. (cs%id_eta_pred_hifreq > 0) .or. (cs%id_etaPF_hifreq > 0) .or. &
2546 (cs%id_uhbt_hifreq > 0) .or. (cs%id_vhbt_hifreq > 0)) &
2547 do_hifreq_output = query_averaging_enabled(cs%diag, time_int_in, time_end_in)
2548 if (do_hifreq_output) then
2549 time_bt_start = time_end_in - real_to_time(dt, unscale=us%T_to_s)
2550 dtbt_diag = dt/(nstep+nfilter) ! Note that this is not dtbt.
2551 endif
2552
2553 ! Zero out the arrays for various time-averaged quantities.
2554 if (find_etaav) then
2555 !$OMP do
2556 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf+1
2557 eta_sum(i,j) = 0.0 ; eta_wtd(i,j) = 0.0
2558 enddo ; enddo
2559 else
2560 !$OMP do
2561 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf+1
2562 eta_wtd(i,j) = 0.0
2563 enddo ; enddo
2564 endif
2565 !$OMP do
2566 do j=js,je ; do i=is-1,ie
2567 cs%ubtav(i,j) = 0.0 ; uhbtav(i,j) = 0.0
2568 pfu_avg(i,j) = 0.0 ; coru_avg(i,j) = 0.0
2569 ldu_avg(i,j) = 0.0 ; ubt_wtd(i,j) = 0.0
2570 enddo ; enddo
2571 !$OMP do
2572 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf
2573 ubt_trans(i,j) = 0.0
2574 enddo ; enddo
2575 !$OMP do
2576 do j=js-1,je ; do i=is,ie
2577 cs%vbtav(i,j) = 0.0 ; vhbtav(i,j) = 0.0
2578 pfv_avg(i,j) = 0.0 ; corv_avg(i,j) = 0.0
2579 ldv_avg(i,j) = 0.0 ; vbt_wtd(i,j) = 0.0
2580 enddo ; enddo
2581 !$OMP do
2582 do j=jsvf-1,jevf ; do i=isvf-1,ievf+1
2583 vbt_trans(i,j) = 0.0
2584 enddo ; enddo
2585 if (integral_bt_cont) then
2586 ubt_int(:,:) = 0.0 ; uhbt_int(:,:) = 0.0
2587 vbt_int(:,:) = 0.0 ; vhbt_int(:,:) = 0.0
2588 endif
2589
2590 p_surf_dyn(:,:) = 0.0
2591 cfl_ltd_vol(:,:) = huge( gv%Z_to_H )
2592 if (cs%bt_limit_integral_transport) then
2593 ! Issue warnings if there are unphysical values of the initial sea surface height or total water column mass.
2594 if (gv%Boussinesq) then
2595 do j=js,je ; do i=is,ie
2596 if ((eta_ic(i,j) < -gv%Z_to_H*g%bathyT(i,j)) .and. (g%mask2dT(i,j) > 0.0)) then
2597 write(mesg,'(ES24.16," vs. ",ES24.16, " at ", ES12.4, ES12.4, i7, i7)') gv%H_to_m*eta(i,j), &
2598 -us%Z_to_m*g%bathyT(i,j), g%geoLonT(i,j), g%geoLatT(i,j), i + g%HI%idg_offset, j + g%HI%jdg_offset
2599 call mom_error(fatal, "btstep: eta_IC starts below bathyT: "//trim(mesg), all_print=.true.)
2600 endif
2601 enddo ; enddo
2602 else
2603 do j=js,je ; do i=is,ie
2604 if ((eta_ic(i,j) < 0.0) .and. (g%mask2dT(i,j) > 0.0)) then
2605 write(mesg,'(" at ", ES12.4, ES12.4, i7, i7)') &
2606 g%geoLonT(i,j), g%geoLatT(i,j), i + g%HI%idg_offset, j + g%HI%jdg_offset
2607 call mom_error(fatal, "btstep: negative eta_IC at start of a non-Boussinesq barotropic solver "//&
2608 trim(mesg), all_print=.true.)
2609 endif
2610 enddo ; enddo
2611 endif
2612 endif
2613
2614 ! Set up the group pass used for halo updates within the barotropic time stepping loops.
2615 call create_group_pass(cs%pass_eta_ubt, eta, cs%BT_Domain)
2616 call create_group_pass(cs%pass_eta_ubt, ubt, vbt, cs%BT_Domain)
2617 if (integral_bt_cont) then
2618 call create_group_pass(cs%pass_eta_ubt, ubt_int, vbt_int, cs%BT_Domain)
2619 ! This is only needed with integral_BT_cont, OBCs and multiple barotropic steps between halo updates.
2620 if (cs%integral_OBCs) &
2621 call create_group_pass(cs%pass_eta_ubt, uhbt_int, vhbt_int, cs%BT_Domain)
2622 endif
2623
2624 ! The following loop contains all of the time steps.
2625 isv = is ; iev = ie ; jsv = js ; jev = je
2626 do n=1,nstep+nfilter
2627 if (cs%clip_velocity) call truncate_velocities(ubt, vbt, dt, g, cs, isv, iev, jsv, jev)
2628
2629 ! Update the range of valid points, either by doing a halo update or by marching inward.
2630 if ((iev - stencil < ie) .or. (jev - stencil < je)) then
2631 if (id_clock_calc > 0) call cpu_clock_end(id_clock_calc)
2632 call do_group_pass(cs%pass_eta_ubt, cs%BT_Domain, clock=id_clock_pass_step)
2633 isv = isvf ; iev = ievf ; jsv = jsvf ; jev = jevf
2634 if (id_clock_calc > 0) call cpu_clock_begin(id_clock_calc)
2635 else
2636 isv = isv+stencil ; iev = iev-stencil
2637 jsv = jsv+stencil ; jev = jev-stencil
2638 endif
2639
2640 ! Store the previous velocities for time-filtered transports and OBCs.
2641 do j=jsv,jev ; do i=isv-2,iev+1
2642 ubt_prev(i,j) = ubt(i,j)
2643 enddo ; enddo
2644 do j=jsv-2,jev+1 ; do i=isv,iev
2645 vbt_prev(i,j) = vbt(i,j)
2646 enddo ; enddo
2647
2648 if (integral_bt_cont) then
2649 !$OMP parallel do default(shared)
2650 do j=jsv-1,jev+1 ; do i=isv-2,iev+1
2651 uhbt_int_prev(i,j) = uhbt_int(i,j)
2652 enddo ; enddo
2653 !$OMP parallel do default(shared)
2654 do j=jsv-2,jev+1 ; do i=isv-1,iev+1
2655 vhbt_int_prev(i,j) = vhbt_int(i,j)
2656 enddo ; enddo
2657 endif
2658
2659 ! Do a predictor step update of eta
2660 if (evolving_face_areas) then
2661 if ((n>1) .and. (mod(n-1,cs%Nonlin_cont_update_period) == 0)) &
2662 call find_face_areas(datu, datv, g, gv, us, cs, ms, 1+iev-ie, eta)
2663 endif
2664
2665 if (cs%dynamic_psurf .or. (.not.cs%BT_project_velocity)) then
2666 ! Estimate the change in the free surface height.
2667 call btloop_eta_predictor(n, dtbt, ubt, vbt, eta, ubt_int, vbt_int, uhbt, vhbt, uhbt0, vhbt0, &
2668 uhbt_int, vhbt_int, btcl_u, btcl_v, datu, datv, eta_ic, eta_src, eta_pred, &
2669 isv, iev, jsv, jev, integral_bt_cont, use_bt_cont, g, us, cs)
2670 endif
2671
2672 if (interp_eta_pf) then
2673 ! Interpolate the effective surface pressure in time
2674 wt_end = n*instep ! This could be (n-0.5)*Instep.
2675 !$OMP parallel do default(shared)
2676 do j=jsv-1,jev+1 ; do i=isv-1,iev+1
2677 eta_pf(i,j) = eta_pf_1(i,j) + wt_end*d_eta_pf(i,j)
2678 enddo ; enddo
2679 endif
2680
2681 v_first = (mod(n+g%first_direction,2)==1)
2682
2683 ! Determine the pressure force accelerations due to the updated eta anomalies.
2684 if (cs%BT_project_velocity) then
2685 call btloop_find_pf(pfu, pfv, isv, iev, jsv, jev, eta, eta_pf, &
2686 gtot_n, gtot_s, gtot_e, gtot_w, dgeo_de, find_etaav, &
2687 wt_accel2(n), eta_sum, v_first, g, us, cs)
2688 else
2689 call btloop_find_pf(pfu, pfv, isv, iev, jsv, jev, eta_pred, eta_pf, &
2690 gtot_n, gtot_s, gtot_e, gtot_w, dgeo_de, find_etaav, &
2691 wt_accel2(n), eta_sum, v_first, g, us, cs)
2692 endif
2693
2694 ! Use the change in eta to determine an additional divergence damping due to the ice strength.
2695 if (cs%dynamic_psurf) then
2696 call btloop_add_dyn_pf(pfu, pfv, eta_pred, eta, dyn_coef_eta, p_surf_dyn, &
2697 isv, iev, jsv, jev, v_first, g, us, cs)
2698 endif
2699
2700 if (v_first) then
2701 ! On odd-steps, update v first.
2702 call btloop_update_v(dtbt, ubt, vbt, v_accel_bt, cor_v, pfv, isv-1, iev+1, jsv-1, jev, &
2703 f_4_v, bt_rem_v, bt_force_v, cor_ref_v, rayleigh_v, &
2704 wt_accel(n), g, us, cs)
2705
2706 ! Now update the zonal velocity.
2707 call btloop_update_u(dtbt, ubt, vbt, u_accel_bt, cor_u, pfu, isv-1, iev, jsv, jev, &
2708 f_4_u, bt_rem_u, bt_force_u, cor_ref_u, rayleigh_u, &
2709 wt_accel(n), g, us, cs)
2710
2711 else
2712 ! On even steps, update u first.
2713 call btloop_update_u(dtbt, ubt, vbt, u_accel_bt, cor_u, pfu, isv-1, iev, jsv-1, jev+1, &
2714 f_4_u, bt_rem_u, bt_force_u, cor_ref_u, rayleigh_u, &
2715 wt_accel(n), g, us, cs)
2716 ! Now update the meridional velocity.
2717 call btloop_update_v(dtbt, ubt, vbt, v_accel_bt, cor_v, pfv, isv, iev, jsv-1, jev, &
2718 f_4_v, bt_rem_v, bt_force_v, cor_ref_v, rayleigh_v, &
2719 wt_accel(n), g, us, cs, cor_bracket_bug=cs%use_old_coriolis_bracket_bug)
2720 endif
2721
2722 ! Determine the transports based on the updated velocities when no OBCs are applied
2723 if (integral_bt_cont) then
2724 if (cs%bt_limit_integral_transport) then
2725 ! Calculate the volume that could be used for divergent transport to use for a limter. This applies to
2726 ! uhbt_int and vhbt_int at each BT step. It does not allow for temporary flooding during the BT cycling.
2727 ! Only the sink is accounted for: if diverent motion occurs at the beginning of the BT cycling but the volume
2728 ! was due only to a +ve source being applied gradually, then the instantaneous eta could drop below the bottom.
2729 if (gv%Boussinesq) then
2730 do j=jsv,jev ; do i=isv,iev
2731 cfl_ltd_vol(i,j) = ( cs%maxCFL_BT_cont * g%areaT(i,j) ) * &
2732 max( 0., ( gv%Z_to_H*g%bathyT(i,j) + eta_ic(i,j) ) + nstep * min( 0., eta_src(i,j) ) )
2733 enddo ; enddo
2734 else
2735 do j=jsv,jev ; do i=isv,iev
2736 cfl_ltd_vol(i,j) = ( cs%maxCFL_BT_cont * g%areaT(i,j) ) * &
2737 max( 0., eta_ic(i,j) + nstep * min( 0., eta_src(i,j) ) )
2738 enddo ; enddo
2739 endif
2740 endif
2741 !$OMP do schedule(static)
2742 do j=jsv,jev ; do i=isv-1,iev
2743 ubt_trans(i,j) = trans_wt1*ubt(i,j) + trans_wt2*ubt_prev(i,j)
2744 ubt_int_prev(i,j) = ubt_int(i,j) ! Store the previous integrated velocity so it can be reset by at OBC points
2745 ubt_int(i,j) = ubt_int(i,j) + dtbt * ubt_trans(i,j)
2746 uhbt_int(i,j) = find_uhbt(ubt_int(i,j), btcl_u(i,j)) + n*dtbt*uhbt0(i,j)
2747 uhbt_int(i,j) = max( -cfl_ltd_vol(i+1,j), min( uhbt_int(i,j), cfl_ltd_vol(i,j) ) )
2748 ! Estimate the mass flux within a single timestep to take the filtered average.
2749 uhbt(i,j) = (uhbt_int(i,j) - uhbt_int_prev(i,j)) * idtbt
2750 enddo ; enddo
2751 !$OMP end do nowait
2752 !$OMP do schedule(static)
2753 do j=jsv-1,jev ; do i=isv,iev
2754 vbt_trans(i,j) = trans_wt1*vbt(i,j) + trans_wt2*vbt_prev(i,j)
2755 vbt_int_prev(i,j) = vbt_int(i,j) ! Store the previous integrated velocity so it can be reset by at OBC points
2756 vbt_int(i,j) = vbt_int(i,j) + dtbt * vbt_trans(i,j)
2757 vhbt_int(i,j) = find_vhbt(vbt_int(i,j), btcl_v(i,j)) + n*dtbt*vhbt0(i,j)
2758 vhbt_int(i,j) = max( -cfl_ltd_vol(i,j+1), min( vhbt_int(i,j), cfl_ltd_vol(i,j) ) )
2759 ! Estimate the mass flux within a single timestep to take the filtered average.
2760 vhbt(i,j) = (vhbt_int(i,j) - vhbt_int_prev(i,j)) * idtbt
2761 enddo ; enddo
2762 elseif (use_bt_cont) then
2763 !$OMP do schedule(static)
2764 do j=jsv,jev ; do i=isv-1,iev
2765 ubt_trans(i,j) = trans_wt1*ubt(i,j) + trans_wt2*ubt_prev(i,j)
2766 uhbt(i,j) = find_uhbt(ubt_trans(i,j), btcl_u(i,j)) + uhbt0(i,j)
2767 enddo ; enddo
2768 !$OMP end do nowait
2769 !$OMP do schedule(static)
2770 do j=jsv-1,jev ; do i=isv,iev
2771 vbt_trans(i,j) = trans_wt1*vbt(i,j) + trans_wt2*vbt_prev(i,j)
2772 vhbt(i,j) = find_vhbt(vbt_trans(i,j), btcl_v(i,j)) + vhbt0(i,j)
2773 enddo ; enddo
2774 else
2775 !$OMP do schedule(static)
2776 do j=jsv,jev ; do i=isv-1,iev
2777 ubt_trans(i,j) = trans_wt1*ubt(i,j) + trans_wt2*ubt_prev(i,j)
2778 uhbt(i,j) = datu(i,j)*ubt_trans(i,j) + uhbt0(i,j)
2779 enddo ; enddo
2780 !$OMP end do nowait
2781 !$OMP do schedule(static)
2782 do j=jsv-1,jev ; do i=isv,iev
2783 vbt_trans(i,j) = trans_wt1*vbt(i,j) + trans_wt2*vbt_prev(i,j)
2784 vhbt(i,j) = datv(i,j)*vbt_trans(i,j) + vhbt0(i,j)
2785 enddo ; enddo
2786 endif
2787
2788 ! This might need to be moved outside of the OMP do loop directives.
2789 if (cs%debug_bt) then
2790 write(mesg,'("BT vel update ",I0)') n
2791 debug_halo = 0 ; if (cs%debug_wide_halos) debug_halo = iev - ie
2792 call uvchksum(trim(mesg)//" PF[uv]", pfu, pfv, cs%debug_BT_HI, haloshift=debug_halo, &
2793 symmetric=.true., unscale=us%L_T_to_m_s*us%s_to_T)
2794 call uvchksum(trim(mesg)//" Cor_[uv]", cor_u, cor_v, cs%debug_BT_HI, haloshift=debug_halo, &
2795 symmetric=.true., unscale=us%L_T_to_m_s*us%s_to_T)
2796 call uvchksum(trim(mesg)//" BT_force_[uv]", bt_force_u, bt_force_v, cs%debug_BT_HI, haloshift=debug_halo, &
2797 symmetric=.true., unscale=us%L_T_to_m_s*us%s_to_T)
2798 call uvchksum(trim(mesg)//" BT_rem_[uv]", bt_rem_u, bt_rem_v, cs%debug_BT_HI, haloshift=debug_halo, &
2799 symmetric=.true., scalar_pair=.true.)
2800 call uvchksum(trim(mesg)//" [uv]bt", ubt, vbt, cs%debug_BT_HI, haloshift=debug_halo, &
2801 symmetric=.true., unscale=us%L_T_to_m_s)
2802 call uvchksum(trim(mesg)//" [uv]bt_trans", ubt_trans, vbt_trans, cs%debug_BT_HI, haloshift=debug_halo, &
2803 symmetric=.true., unscale=us%L_T_to_m_s)
2804 call uvchksum(trim(mesg)//" [uv]hbt", uhbt, vhbt, cs%debug_BT_HI, haloshift=debug_halo, &
2805 symmetric=.true., unscale=us%s_to_T*us%L_to_m**2*gv%H_to_m)
2806 if (integral_bt_cont) &
2807 call uvchksum(trim(mesg)//" [uv]hbt_int", uhbt_int, vhbt_int, cs%debug_BT_HI, haloshift=debug_halo, &
2808 symmetric=.true., unscale=us%L_to_m**2*gv%H_to_m)
2809 endif
2810
2811 ! Apply open boundary condition considerations to revise the updated velocities and transports.
2812 if (cs%BT_OBC%u_OBCs_on_PE) then
2813 !$OMP single
2814 call apply_u_velocity_obcs(ubt, uhbt, ubt_trans, eta, spv_col_avg, ubt_prev, bt_obc, &
2815 g, ms, gv, us, cs, iev-ie, dtbt, cs%bebt, use_bt_cont, integral_bt_cont, n*dtbt, &
2816 datu, btcl_u, uhbt0, ubt_int, ubt_int_prev, uhbt_int, uhbt_int_prev)
2817 !$OMP end single
2818 endif
2819
2820 if (cs%BT_OBC%v_OBCs_on_PE) then
2821 !$OMP single
2822 call apply_v_velocity_obcs(vbt, vhbt, vbt_trans, eta, spv_col_avg, vbt_prev, bt_obc, &
2823 g, ms, gv, us, cs, iev-ie, dtbt, cs%bebt, use_bt_cont, integral_bt_cont, n*dtbt, &
2824 datv, btcl_v, vhbt0, vbt_int, vbt_int_prev, vhbt_int, vhbt_int_prev)
2825 !$OMP end single
2826 endif
2827
2828 ! Contribute to the running sums of the transports and velocities.
2829 !$OMP do
2830 do j=js,je ; do i=is-1,ie
2831 cs%ubtav(i,j) = cs%ubtav(i,j) + wt_trans(n) * ubt_trans(i,j)
2832 uhbtav(i,j) = uhbtav(i,j) + wt_trans(n) * uhbt(i,j)
2833 ubt_wtd(i,j) = ubt_wtd(i,j) + wt_vel(n) * ubt(i,j)
2834 enddo ; enddo
2835 !$OMP end do nowait
2836 !$OMP do
2837 do j=js-1,je ; do i=is,ie
2838 cs%vbtav(i,j) = cs%vbtav(i,j) + wt_trans(n) * vbt_trans(i,j)
2839 vhbtav(i,j) = vhbtav(i,j) + wt_trans(n) * vhbt(i,j)
2840 vbt_wtd(i,j) = vbt_wtd(i,j) + wt_vel(n) * vbt(i,j)
2841 enddo ; enddo
2842 !$OMP end do nowait
2843
2844 if (cs%debug_bt) then
2845 call uvchksum("BT [uv]hbt just after OBC", uhbt, vhbt, cs%debug_BT_HI, haloshift=debug_halo, &
2846 symmetric=.true., unscale=us%s_to_T*us%L_to_m**2*gv%H_to_m)
2847 if (integral_bt_cont) &
2848 call uvchksum("BT [uv]hbt_int just after OBC", uhbt_int, vhbt_int, cs%debug_BT_HI, &
2849 haloshift=debug_halo, symmetric=.true., unscale=us%L_to_m**2*gv%H_to_m)
2850 endif
2851
2852 ! Update eta in a corrector step using the barotropic continuity equation.
2853 if (integral_bt_cont) then
2854 eta_cor_multiplier = n
2855 if ( cs%bt_adjust_src_for_filter ) then
2856 if ( nstep > nfilter ) then
2857 eta_cor_multiplier = min(nstep - nfilter, n) * nstep / real(nstep - nfilter)
2858 else
2859 eta_cor_multiplier = nstep
2860 endif
2861 endif
2862 !$OMP do
2863 do j=jsv,jev ; do i=isv,iev
2864 eta(i,j) = (eta_ic(i,j) + eta_cor_multiplier*eta_src(i,j)) + cs%IareaT_OBCmask(i,j) * &
2865 ((uhbt_int(i-1,j) - uhbt_int(i,j)) + (vhbt_int(i,j-1) - vhbt_int(i,j)))
2866 ! eta_acc contains the magnitude of the largest term in the above expression which
2867 ! will be used to estimate a bound for round off when comparing to the bottom depth
2868 eta_acc = abs( cs%IareaT_OBCmask(i,j) * &
2869 ((uhbt_int(i-1,j) - uhbt_int(i,j)) + (vhbt_int(i,j-1) - vhbt_int(i,j))) )
2870 eta_acc = max( eta_acc, abs( eta_cor_multiplier*eta_src(i,j) ), abs( eta_ic(i,j) ) )
2871 if ( g%mask2dT(i,j) * ( eta(i,j) + gv%Z_to_H*g%bathyT(i,j) ) > &
2872 -g%mask2dT(i,j) * eta_acc * epsilon(eta_acc) * 2. ) &
2873 eta(i,j) = max( eta(i,j), -gv%Z_to_H*g%bathyT(i,j) )
2874 eta_wtd(i,j) = eta_wtd(i,j) + eta(i,j) * wt_eta(n)
2875 if ((eta(i,j) < -gv%Z_to_H*g%bathyT(i,j)) .and. (g%mask2dT(i,j) > 0.0)) then
2876 write(mesg,'(ES24.16," vs. ",ES24.16, " at ", ES12.4, ES12.4, i7, i7)') gv%H_to_m*eta(i,j), &
2877 -us%Z_to_m*g%bathyT(i,j), g%geoLonT(i,j), g%geoLatT(i,j), i + g%HI%idg_offset, j + g%HI%jdg_offset
2878 if (cs%bt_limit_integral_transport) &
2879 call mom_error(fatal, "btstep: eta has dropped below bathyT: "//trim(mesg))
2880 endif
2881 enddo ; enddo
2882 else
2883 !$OMP do
2884 do j=jsv,jev ; do i=isv,iev
2885 eta(i,j) = (eta(i,j) + eta_src(i,j)) + (dtbt * cs%IareaT_OBCmask(i,j)) * &
2886 ((uhbt(i-1,j) - uhbt(i,j)) + (vhbt(i,j-1) - vhbt(i,j)))
2887 eta_wtd(i,j) = eta_wtd(i,j) + eta(i,j) * wt_eta(n)
2888 enddo ; enddo
2889 endif
2890
2891 if (cs%debug_bt) then
2892 write(mesg,'("BT step ",I0)') n
2893 call uvchksum(trim(mesg)//" [uv]bt", ubt, vbt, cs%debug_BT_HI, haloshift=debug_halo, &
2894 symmetric=.true., unscale=us%L_T_to_m_s)
2895 call hchksum(eta, trim(mesg)//" eta", cs%debug_BT_HI, haloshift=debug_halo, unscale=gv%H_to_MKS)
2896 endif
2897
2898 ! Issue warnings if there are unphysical values of the sea surface height or total water column mass.
2899 if (gv%Boussinesq) then
2900 do j=js,je ; do i=is,ie
2901 if ((eta(i,j) < -gv%Z_to_H*g%bathyT(i,j)) .and. (g%mask2dT(i,j) > 0.0)) then
2902 write(mesg,'(ES24.16," vs. ",ES24.16, " at ", ES12.4, ES12.4, i7, i7)') gv%H_to_m*eta(i,j), &
2903 -us%Z_to_m*g%bathyT(i,j), g%geoLonT(i,j), g%geoLatT(i,j), i + g%HI%idg_offset, j + g%HI%jdg_offset
2904 if (cs%bt_limit_integral_transport) &
2905 call mom_error(fatal, "btstep: eta has dropped below bathyT: "//trim(mesg))
2906 if (err_count < 2) &
2907 call mom_error(warning, "btstep: eta has dropped below bathyT: "//trim(mesg), all_print=.true.)
2908 err_count = err_count + 1
2909 endif
2910 enddo ; enddo
2911 else
2912 do j=js,je ; do i=is,ie
2913 if ((eta(i,j) < 0.0) .and. (g%mask2dT(i,j) > 0.0)) then
2914 write(mesg,'(" at ", ES12.4, ES12.4, i7, i7)') &
2915 g%geoLonT(i,j), g%geoLatT(i,j), i + g%HI%idg_offset, j + g%HI%jdg_offset
2916 if (cs%bt_limit_integral_transport) &
2917 call mom_error(fatal, "btstep: negative eta in a non-Boussinesq barotropic solver "//trim(mesg))
2918 if (err_count < 2) &
2919 call mom_error(warning, "btstep: negative eta in a non-Boussinesq barotropic solver "//&
2920 trim(mesg), all_print=.true.)
2921 err_count = err_count + 1
2922 endif
2923 enddo ; enddo
2924 endif
2925
2926 ! Accumulate some diagnostics of time-averaged barotropic accelerations.
2927 if (do_ave) then
2928 if ((cs%id_PFu_bt > 0) .or. associated(adp%bt_pgf_u)) then
2929 !$OMP do
2930 do j=js,je ; do i=is-1,ie
2931 pfu_avg(i,j) = pfu_avg(i,j) + wt_accel2(n) * pfu(i,j)
2932 enddo ; enddo
2933 !$OMP end do nowait
2934 endif
2935 if ((cs%id_PFv_bt > 0) .or. associated(adp%bt_pgf_v)) then
2936 !$OMP do
2937 do j=js-1,je ; do i=is,ie
2938 pfv_avg(i,j) = pfv_avg(i,j) + wt_accel2(n) * pfv(i,j)
2939 enddo ; enddo
2940 !$OMP end do nowait
2941 endif
2942 if ((cs%id_Coru_bt > 0) .or. associated(adp%bt_cor_u)) then
2943 !$OMP do
2944 do j=js,je ; do i=is-1,ie
2945 coru_avg(i,j) = coru_avg(i,j) + wt_accel2(n) * cor_u(i,j)
2946 enddo ; enddo
2947 !$OMP end do nowait
2948 endif
2949 if ((cs%id_Corv_bt > 0) .or. associated(adp%bt_cor_v)) then
2950 !$OMP do
2951 do j=js-1,je ; do i=is,ie
2952 corv_avg(i,j) = corv_avg(i,j) + wt_accel2(n) * cor_v(i,j)
2953 enddo ; enddo
2954 !$OMP end do nowait
2955 endif
2956
2957 if (cs%linear_wave_drag) then
2958 if ((cs%id_LDu_bt > 0) .or. (associated(adp%bt_lwd_u))) then
2959 !$OMP do
2960 do j=js,je ; do i=is-1,ie
2961 ldu_avg(i,j) = ldu_avg(i,j) - wt_accel2(n) * (ubt(i,j) * rayleigh_u(i,j))
2962 enddo ; enddo
2963 !$OMP end do nowait
2964 endif
2965 if ((cs%id_LDv_bt > 0) .or. (associated(adp%bt_lwd_v))) then
2966 !$OMP do
2967 do j=js-1,je ; do i=is,ie
2968 ldv_avg(i,j) = ldv_avg(i,j) - wt_accel2(n) * (vbt(i,j) * rayleigh_v(i,j))
2969 enddo ; enddo
2970 !$OMP end do nowait
2971 endif
2972 endif
2973 endif
2974
2975 if (do_hifreq_output) then
2976 ! Note that this compresses the time so that all of the timesteps, including those in the
2977 ! extra timesteps for filtering, fit within dt.
2978 time_step_end = time_bt_start + real_to_time(n*dtbt_diag, unscale=us%T_to_s)
2979 call enable_averages(dtbt, time_step_end, cs%diag)
2980 if (cs%id_ubt_hifreq > 0) call post_data(cs%id_ubt_hifreq, ubt(isdb:iedb,jsd:jed), cs%diag)
2981 if (cs%id_vbt_hifreq > 0) call post_data(cs%id_vbt_hifreq, vbt(isd:ied,jsdb:jedb), cs%diag)
2982 if (cs%id_eta_hifreq > 0) call post_data(cs%id_eta_hifreq, eta(isd:ied,jsd:jed), cs%diag)
2983 if (cs%id_etaPF_hifreq > 0) then
2984 if (cs%BT_project_velocity) then
2985 do j=js,je ; do i=is,ie
2986 eta_anom_pf(i,j) = eta(i,j) - eta_pf(i,j)
2987 enddo ; enddo
2988 else
2989 do j=js,je ; do i=is,ie
2990 eta_anom_pf(i,j) = eta_pred(i,j) - eta_pf(i,j)
2991 enddo ; enddo
2992 endif
2993 call post_data(cs%id_etaPF_hifreq, eta_anom_pf(isd:ied,jsd:jed), cs%diag)
2994 endif
2995 if (cs%id_uhbt_hifreq > 0) call post_data(cs%id_uhbt_hifreq, uhbt(isdb:iedb,jsd:jed), cs%diag)
2996 if (cs%id_vhbt_hifreq > 0) call post_data(cs%id_vhbt_hifreq, vhbt(isd:ied,jsdb:jedb), cs%diag)
2997 if (cs%BT_project_velocity) then
2998 ! This diagnostic is redundant in this case and should probably be omitted.
2999 if (cs%id_eta_pred_hifreq > 0) call post_data(cs%id_eta_pred_hifreq, eta(isd:ied,jsd:jed), cs%diag)
3000 else
3001 if (cs%id_eta_pred_hifreq > 0) call post_data(cs%id_eta_pred_hifreq, eta_pred(isd:ied,jsd:jed), cs%diag)
3002 endif
3003 endif
3004 enddo ! end of do n=1,ntimestep
3005
3006 ! Reset the time information in the diag type.
3007 if (do_hifreq_output) call enable_averaging(time_int_in, time_end_in, cs%diag)
3008
3009end subroutine btstep_timeloop
3010
3011
3012!> Find the Coriolis force terms _zon and _mer.
3013subroutine btstep_find_cor(q, DCor_u, DCor_v, f_4_u, f_4_v, isvf, ievf, jsvf, jevf, CS)
3014 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3015 real, intent(in) :: q(SZIBW_(CS),SZJBW_(CS)) !< A pseudo potential vorticity [T-1 Z-1 ~> s-1 m-1]
3016 !! or [T-1 H-1 ~> s-1 m-1 or m2 s-1 kg-1]
3017 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3018 DCor_u !< An averaged depth or total thickness at u points [Z ~> m] or [H ~> m or kg m-2].
3019 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3020 DCor_v !< An averaged depth or total thickness at v points [Z ~> m] or [H ~> m or kg m-2].
3021 real, dimension(4,SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3022 f_4_u !< The terms giving the contribution to the Coriolis acceleration at a zonal
3023 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
3024 !! These are the products of thicknesses at v points and appropriately staggered
3025 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
3026 !! they are approximately f / 4. The 4 values on the innermost loop are for
3027 !! v-velocities to the southwest, southeast, northwest and northeast.
3028 real, dimension(4,SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3029 f_4_v !< The terms giving the contribution to the Coriolis acceleration at a meridional
3030 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
3031 !! These are the products of thicknesses at u points and appropriately staggered
3032 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
3033 !! they are approximately f / 4. The 4 values on the innermost loop are for
3034 !! u-velocities to the southwest, southeast, northwest and northeast.
3035 integer, intent(in) :: isvf !< The starting i-index of the largest valid range for tracer points
3036 integer, intent(in) :: ievf !< The ending i-index of the largest valid range for tracer points
3037 integer, intent(in) :: jsvf !< The starting j-index of the largest valid range for tracer points
3038 integer, intent(in) :: jevf !< The ending j-index of the largest valid range for tracer points
3039
3040 ! real :: C1_3 ! One third [nondim]
3041 integer :: i, j
3042
3043 if (cs%Sadourny) then
3044 !$OMP parallel do default(shared)
3045 do j=jsvf-1,jevf ; do i=isvf-1,ievf+1
3046 f_4_v(1,i,j) = cs%OBCmask_v(i,j) * dcor_u(i-1,j) * q(i-1,j)
3047 f_4_v(2,i,j) = cs%OBCmask_v(i,j) * dcor_u(i,j) * q(i,j)
3048 f_4_v(4,i,j) = cs%OBCmask_v(i,j) * dcor_u(i,j+1) * q(i,j)
3049 f_4_v(3,i,j) = cs%OBCmask_v(i,j) * dcor_u(i-1,j+1) * q(i-1,j)
3050 enddo ; enddo
3051 !$OMP parallel do default(shared)
3052 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf
3053 f_4_u(4,i,j) = cs%OBCmask_u(i,j) * dcor_v(i+1,j) * q(i,j)
3054 f_4_u(3,i,j) = cs%OBCmask_u(i,j) * dcor_v(i,j) * q(i,j)
3055 f_4_u(1,i,j) = cs%OBCmask_u(i,j) * dcor_v(i,j-1) * q(i,j-1)
3056 f_4_u(2,i,j) = cs%OBCmask_u(i,j) * dcor_v(i+1,j-1) * q(i,j-1)
3057 enddo ; enddo
3058 else !### if (CS%answer_date < 20250601) then ! Uncomment this later.
3059 !$OMP parallel do default(shared)
3060 do j=jsvf-1,jevf ; do i=isvf-1,ievf+1
3061 f_4_v(1,i,j) = cs%OBCmask_v(i,j) * dcor_u(i-1,j) * ((q(i,j) + q(i-1,j-1)) + q(i-1,j)) / 3.0
3062 f_4_v(2,i,j) = cs%OBCmask_v(i,j) * dcor_u(i,j) * (q(i,j) + (q(i-1,j) + q(i,j-1))) / 3.0
3063 f_4_v(4,i,j) = cs%OBCmask_v(i,j) * dcor_u(i,j+1) * (q(i,j) + (q(i-1,j) + q(i,j+1))) / 3.0
3064 f_4_v(3,i,j) = cs%OBCmask_v(i,j) * dcor_u(i-1,j+1) * ((q(i,j) + q(i-1,j+1)) + q(i-1,j)) / 3.0
3065 enddo ; enddo
3066 !$OMP parallel do default(shared)
3067 do j=jsvf-1,jevf+1 ; do i=isvf-1,ievf
3068 f_4_u(4,i,j) = cs%OBCmask_u(i,j) * dcor_v(i+1,j) * (q(i,j) + (q(i+1,j) + q(i,j-1))) / 3.0
3069 f_4_u(3,i,j) = cs%OBCmask_u(i,j) * dcor_v(i,j) * (q(i,j) + (q(i-1,j) + q(i,j-1))) / 3.0
3070 f_4_u(1,i,j) = cs%OBCmask_u(i,j) * dcor_v(i,j-1) * ((q(i,j) + q(i-1,j-1)) + q(i,j-1)) / 3.0
3071 f_4_u(2,i,j) = cs%OBCmask_u(i,j) * dcor_v(i+1,j-1) * ((q(i,j) + q(i+1,j-1)) + q(i,j-1)) / 3.0
3072 enddo ; enddo
3073 ! else
3074 ! C1_3 = 1.0 / 3.0
3075 ! !$OMP parallel do default(shared)
3076 ! do J=jsvf-1,jevf ; do i=isvf-1,ievf+1
3077 ! f_4_v(1,i,J) = CS%OBCmask_v(i,J) * DCor_u(I-1,j) * ((q(I,J) + q(I-1,J-1)) + q(I-1,J)) * C1_3
3078 ! f_4_v(2,i,J) = CS%OBCmask_v(i,J) * DCor_u(I,j) * (q(I,J) + (q(I-1,J) + q(I,J-1))) * C1_3
3079 ! f_4_v(4,i,J) = CS%OBCmask_v(i,J) * DCor_u(I,j+1) * (q(I,J) + (q(I-1,J) + q(I,J+1))) * C1_3
3080 ! f_4_v(3,i,J) = CS%OBCmask_v(i,J) * DCor_u(I-1,j+1) * ((q(I,J) + q(I-1,J+1)) + q(I-1,J)) * C1_3
3081 ! enddo ; enddo
3082 ! !$OMP parallel do default(shared)
3083 ! do j=jsvf-1,jevf+1 ; do I=isvf-1,ievf
3084 ! f_4_u(4,I,j) = CS%OBCmask_u(I,j) * DCor_v(i+1,J) * (q(I,J) + (q(I+1,J) + q(I,J-1))) * C1_3
3085 ! f_4_u(3,I,j) = CS%OBCmask_u(I,j) * DCor_v(i,J) * (q(I,J) + (q(I-1,J) + q(I,J-1))) * C1_3
3086 ! f_4_u(1,I,j) = CS%OBCmask_u(I,j) * DCor_v(i,J-1) * ((q(I,J) + q(I-1,J-1)) + q(I,J-1)) * C1_3
3087 ! f_4_u(2,I,j) = CS%OBCmask_u(I,j) * DCor_v(i+1,J-1) * ((q(I,J) + q(I+1,J-1)) + q(I,J-1)) * C1_3
3088 ! enddo ; enddo
3089 endif
3090
3091end subroutine btstep_find_cor
3092
3093!> Do a CFL-based truncation of any excessively large batotropic velocities.
3094!! This should only be used as desperate debugging measure.
3095subroutine truncate_velocities(ubt, vbt, dt, G, CS, isv, iev, jsv, jev)
3096 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3097 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3098 real, intent(inout) :: ubt(SZIBW_(CS),SZJW_(CS)) !< The zonal barotropic velocity [L T-1 ~> m s-1]
3099 real, intent(inout) :: vbt(SZIW_(CS),SZJBW_(CS)) !< The meridional barotropic velocity [L T-1 ~> m s-1]
3100 real, intent(in) :: dt !< The time increment to integrate over [T ~> s].
3101 integer, intent(in) :: isv !< The starting valid tracer array i-index that is being worked on
3102 integer, intent(in) :: iev !< The ending valid tracer array i-index that is being worked on
3103 integer, intent(in) :: jsv !< The starting valid tracer array j-index that is being worked on
3104 integer, intent(in) :: jev !< The ending valid tracer array j-index being that is worked on
3105
3106 integer :: i, j
3107
3108 if (cs%clip_velocity) then
3109 do j=jsv,jev ; do i=isv-1,iev
3110 if ((ubt(i,j) * (dt * g%dy_Cu(i,j))) * g%IareaT(i+1,j) < -cs%CFL_trunc) then
3111 ! Add some error reporting later.
3112 ubt(i,j) = (-0.95*cs%CFL_trunc) * (g%areaT(i+1,j) / (dt * g%dy_Cu(i,j)))
3113 elseif ((ubt(i,j) * (dt * g%dy_Cu(i,j))) * g%IareaT(i,j) > cs%CFL_trunc) then
3114 ! Add some error reporting later.
3115 ubt(i,j) = (0.95*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dy_Cu(i,j)))
3116 endif
3117 enddo ; enddo
3118 do j=jsv-1,jev ; do i=isv,iev
3119 if ((vbt(i,j) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j+1) < -cs%CFL_trunc) then
3120 ! Add some error reporting later.
3121 vbt(i,j) = (-0.9*cs%CFL_trunc) * (g%areaT(i,j+1) / (dt * g%dx_Cv(i,j)))
3122 elseif ((vbt(i,j) * (dt * g%dx_Cv(i,j))) * g%IareaT(i,j) > cs%CFL_trunc) then
3123 ! Add some error reporting later.
3124 vbt(i,j) = (0.9*cs%CFL_trunc) * (g%areaT(i,j) / (dt * g%dx_Cv(i,j)))
3125 endif
3126 enddo ; enddo
3127 endif
3128
3129end subroutine truncate_velocities
3130
3131
3132!> A routine to set eta_pred and the running time integral of uhbt and vhbt.
3133subroutine btloop_eta_predictor(n, dtbt, ubt, vbt, eta, ubt_int, vbt_int, uhbt, vhbt, uhbt0, vhbt0, &
3134 uhbt_int, vhbt_int, BTCL_u, BTCL_v, Datu, Datv, &
3135 eta_IC, eta_src, eta_pred, isv, iev, jsv, jev, &
3136 integral_BT_cont, use_BT_cont, G, US, CS)
3137 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
3138 type(barotropic_cs), intent(in) :: CS !< Barotropic control structure
3139 integer, intent(in) :: n !< The current step in loop of timesteps
3140 real, intent(in) :: dtbt !< The barotropic time step [T ~> s]
3141 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3142 ubt !< The zonal barotropic velocity [L T-1 ~> m s-1].
3143 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3144 vbt !< The zonal barotropic velocity [L T-1 ~> m s-1].
3145 real, target, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3146 eta !< The barotropic free surface height anomaly or column mass
3147 !! anomaly [H ~> m or kg m-2]
3148 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3149 ubt_int !< The running time integral of ubt over the time steps [L ~> m].
3150 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3151 vbt_int !< The running time integral of vbt over the time steps [L ~> m].
3152 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3153 uhbt0 !< The difference between the sum of the layer zonal thickness
3154 !! fluxes and the barotropic thickness flux using the same
3155 !! velocity [H L2 T-1 ~> m3 s-1 or kg s-1].
3156 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3157 vhbt0 !< The difference between the sum of the layer meridional
3158 !! thickness fluxes and the barotropic thickness flux using
3159 !! the same velocities [H L2 T-1 ~> m3 s-1 or kg s-1].
3160 type(local_bt_cont_u_type), dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3161 BTCL_u !< A repackaged version of the u-point information in BT_cont.
3162 type(local_bt_cont_v_type), dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3163 BTCL_v !< A repackaged version of the v-point information in BT_cont.
3164 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3165 Datu !< Basin depth at u-velocity grid points times the y-grid
3166 !! spacing [H L ~> m2 or kg m-1].
3167 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3168 Datv !< Basin depth at v-velocity grid points times the x-grid
3169 !! spacing [H L ~> m2 or kg m-1].
3170 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3171 eta_IC !< A local copy of the initial 2-D eta field (eta_in) [H ~> m or kg m-2]
3172 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3173 eta_src !< The source of eta per barotropic timestep [H ~> m or kg m-2].
3174 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3175 uhbt !< The zonal barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].
3176 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3177 vhbt !< The meridional barotropic thickness fluxes [H L2 T-1 ~> m3 s-1 or kg s-1].
3178 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3179 uhbt_int !< The running time integral of uhbt over the time steps [H L2 ~> m3 or kg].
3180 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3181 vhbt_int !< The running time integral of vhbt over the time steps [H L2 ~> m3 or kg].
3182 real, target, dimension(SZIW_(CS),SZJW_(CS)), intent(inout) :: &
3183 eta_pred !< A predictor value of eta [H ~> m or kg m-2] like eta.
3184 integer, intent(in) :: isv !< The starting i-index of eta_pred to calculate
3185 integer, intent(in) :: iev !< The ending i-index of eta_pred to calculate
3186 integer, intent(in) :: jsv !< The starting j-index of eta_pred to calculate
3187 integer, intent(in) :: jev !< The ending j-index of eta_pred to calculate
3188 logical, intent(in) :: integral_BT_cont !< If true, update the barotropic continuity equation directly
3189 !! from the initial condition using the time-integrated barotropic velocity.
3190 logical, intent(in) :: use_BT_cont !< If true, use the information in the BT_cont_type to determine
3191 !! barotropic transports as a function of the barotropic velocities.
3192 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3193
3194 integer :: i, j
3195
3196 !$OMP parallel default(shared)
3197 if (integral_bt_cont) then
3198 !$OMP do
3199 do j=jsv-1,jev+1 ; do i=isv-2,iev+1
3200 uhbt_int(i,j) = find_uhbt(ubt_int(i,j) + dtbt*ubt(i,j), btcl_u(i,j)) + n*dtbt*uhbt0(i,j)
3201 enddo ; enddo
3202 !$OMP end do nowait
3203 !$OMP do
3204 do j=jsv-2,jev+1 ; do i=isv-1,iev+1
3205 vhbt_int(i,j) = find_vhbt(vbt_int(i,j) + dtbt*vbt(i,j), btcl_v(i,j)) + n*dtbt*vhbt0(i,j)
3206 enddo ; enddo
3207 !$OMP do
3208 do j=jsv-1,jev+1 ; do i=isv-1,iev+1
3209 eta_pred(i,j) = (eta_ic(i,j) + n*eta_src(i,j)) + cs%IareaT_OBCmask(i,j) * &
3210 ((uhbt_int(i-1,j) - uhbt_int(i,j)) + (vhbt_int(i,j-1) - vhbt_int(i,j)))
3211 enddo ; enddo
3212 elseif (use_bt_cont) then
3213 !$OMP do
3214 do j=jsv-1,jev+1 ; do i=isv-2,iev+1
3215 uhbt(i,j) = find_uhbt(ubt(i,j), btcl_u(i,j)) + uhbt0(i,j)
3216 enddo ; enddo
3217 !$OMP do
3218 do j=jsv-2,jev+1 ; do i=isv-1,iev+1
3219 vhbt(i,j) = find_vhbt(vbt(i,j), btcl_v(i,j)) + vhbt0(i,j)
3220 enddo ; enddo
3221 !$OMP do
3222 do j=jsv-1,jev+1 ; do i=isv-1,iev+1
3223 eta_pred(i,j) = (eta(i,j) + eta_src(i,j)) + (dtbt * cs%IareaT_OBCmask(i,j)) * &
3224 ((uhbt(i-1,j) - uhbt(i,j)) + (vhbt(i,j-1) - vhbt(i,j)))
3225 enddo ; enddo
3226 else
3227 !$OMP do
3228 do j=jsv-1,jev+1 ; do i=isv-1,iev+1
3229 eta_pred(i,j) = (eta(i,j) + eta_src(i,j)) + (dtbt * cs%IareaT_OBCmask(i,j)) * &
3230 (((datu(i-1,j)*ubt(i-1,j) + uhbt0(i-1,j)) - &
3231 (datu(i,j)*ubt(i,j) + uhbt0(i,j))) + &
3232 ((datv(i,j-1)*vbt(i,j-1) + vhbt0(i,j-1)) - &
3233 (datv(i,j)*vbt(i,j) + vhbt0(i,j))))
3234 enddo ; enddo
3235 endif
3236 !$OMP end parallel
3237
3238end subroutine btloop_eta_predictor
3239
3240subroutine btloop_find_pf(PFu, PFv, isv, iev, jsv, jev, eta_PF_BT, eta_PF, &
3241 gtot_N, gtot_S, gtot_E, gtot_W, dgeo_de, find_etaav, &
3242 wt_accel2_n, eta_sum, v_first, G, US, CS)
3243 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3244 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3245 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3246 PFu !< The anomalous zonal pressure force acceleration [L T-2 ~> m s-2].
3247 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3248 PFv !< The meridional pressure force acceleration [L T-2 ~> m s-2].
3249 integer, intent(in) :: isv !< The starting i-index of eta being set in ths loop
3250 integer, intent(in) :: iev !< The ending i-index of eta_pred being set in ths loop
3251 integer, intent(in) :: jsv !< The starting j-index of eta_pred being set in ths loop
3252 integer, intent(in) :: jev !< The ending j-index of eta_pred being set in ths loop
3253 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3254 eta_PF_BT !< The eta array (either the SSH anomaly or column mass anomaly) that
3255 !! determines the barotropic pressure force [H ~> m or kg m-2]
3256 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3257 eta_PF !< The input 2-D eta field (either SSH anomaly or column mass anomaly)
3258 !! that was used to calculate the input pressure gradient
3259 !! accelerations [H ~> m or kg m-2].
3260 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3261 gtot_N !< The effective total reduced gravity used to relate free surface height
3262 !! deviations to pressure forces (including GFS and baroclinic contributions)
3263 !! in the barotropic momentum equations half a grid-point to the north of a
3264 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3265 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3266 gtot_S !< The effective total reduced gravity used to relate free surface height
3267 !! deviations to pressure forces (including GFS and baroclinic contributions)
3268 !! in the barotropic momentum equations half a grid-point to the south of a
3269 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3270 !! (See Hallberg, J Comp Phys 1997 for a discussion of gtot_E and gtot_W.)
3271 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3272 gtot_E !< The effective total reduced gravity used to relate free surface height
3273 !! deviations to pressure forces (including GFS and baroclinic contributions)
3274 !! in the barotropic momentum equations half a grid-point to the east of a
3275 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3276 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3277 gtot_W !< The effective total reduced gravity used to relate free surface height
3278 !! deviations to pressure forces (including GFS and baroclinic contributions)
3279 !! in the barotropic momentum equations half a grid-point to the west of a
3280 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3281 !! (See Hallberg, J Comp Phys 1997 for a discussion of gtot_E and gtot_W.)
3282 real, intent(in) :: dgeo_de !< The constant of proportionality between geopotential and
3283 !! sea surface height [nondim]. It is of order 1, but for stability this
3284 !! may be made larger than the physical problem would suggest.
3285 logical, intent(in) :: find_etaav !< If true, diagnose the time mean value of eta
3286 real, intent(in) :: wt_accel2_n !< The weighting value of wt_accel2 at step n.
3287 real, dimension(SZIW_(CS),SZJW_(CS)), intent(inout) :: &
3288 eta_sum !< A weighted running sum of eta summed across the timesteps [H ~> m or kg m-2]
3289 logical, intent(in) :: v_first !< If true, update the v-velocity first with the present loop iteration
3290 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3291
3292 ! Local variables
3293 integer :: i, j, js_u, je_u, is_v, ie_v
3294
3295 ! Ensure that the extra points used for the temporally staggered Coriolis terms are updated.
3296 if (v_first) then
3297 is_v = isv-1 ; ie_v = iev+1 ; js_u = jsv ; je_u = jev
3298 else
3299 is_v = isv ; ie_v = iev ; js_u = jsv-1 ; je_u = jev+1
3300 endif
3301
3302 !$OMP do schedule(static)
3303 do j=js_u,je_u ; do i=isv-1,iev
3304 pfu(i,j) = (((eta_pf_bt(i,j)-eta_pf(i,j))*gtot_e(i,j)) - &
3305 ((eta_pf_bt(i+1,j)-eta_pf(i+1,j))*gtot_w(i+1,j))) * &
3306 dgeo_de * cs%IdxCu(i,j)
3307 enddo ; enddo
3308 !$OMP end do nowait
3309
3310 !$OMP do schedule(static)
3311 do j=jsv-1,jev ; do i=is_v,ie_v
3312 pfv(i,j) = (((eta_pf_bt(i,j)-eta_pf(i,j))*gtot_n(i,j)) - &
3313 ((eta_pf_bt(i,j+1)-eta_pf(i,j+1))*gtot_s(i,j+1))) * &
3314 dgeo_de * cs%IdyCv(i,j)
3315 enddo ; enddo
3316 !$OMP end do nowait
3317
3318 if (find_etaav .and. (abs(wt_accel2_n) > 0.0)) then
3319 !$OMP do
3320 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3321 eta_sum(i,j) = eta_sum(i,j) + wt_accel2_n * eta_pf_bt(i,j)
3322 enddo ; enddo
3323 !$OMP end do nowait
3324 endif
3325
3326end subroutine btloop_find_pf
3327
3328!> This routine adds a dynamic pressure force based on the temporal changes in the predicted value
3329!! of eta, perhaps as an effective divergence damping to emulate the rigidity of an ice-sheet.
3330subroutine btloop_add_dyn_pf(PFu, PFv, eta_pred, eta, dyn_coef_eta, p_surf_dyn, &
3331 isv, iev, jsv, jev, v_first, G, US, CS)
3332 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3333 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3334 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3335 PFu !< The anomalous zonal pressure force acceleration [L T-2 ~> m s-2].
3336 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3337 PFv !< The meridional pressure force acceleration [L T-2 ~> m s-2].
3338 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3339 eta_pred !< The updated eta field (either SSH anomaly or column mass anomaly) that is
3340 !! used to estimate the divergence that is to be damped [H ~> m or kg m-2].
3341 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3342 eta !< The previous eta field (either SSH anomaly or column mass anomaly) that is
3343 !! used to estimate the divergence that is to be damped [H ~> m or kg m-2].
3344 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3345 dyn_coef_eta !< The coefficient relating the changes in eta to the dynamic surface pressure
3346 !! under rigid ice [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
3347 real, dimension(SZIW_(CS),SZJW_(CS)), intent(inout) :: &
3348 p_surf_dyn !< A dynamic surface pressure under rigid ice [L2 T-2 ~> m2 s-2].
3349 integer, intent(in) :: isv !< The starting i-index of eta being set in ths loop
3350 integer, intent(in) :: iev !< The ending i-index of eta_pred being set in ths loop
3351 integer, intent(in) :: jsv !< The starting j-index of eta_pred being set in ths loop
3352 integer, intent(in) :: jev !< The ending j-index of eta_pred being set in ths loop
3353 logical, intent(in) :: v_first !< If true, update the v-velocity first with the present loop iteration
3354 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3355
3356 ! Local variables
3357 integer :: i, j, js_u, je_u, is_v, ie_v
3358
3359 ! Ensure that the extra points used for the temporally staggered Coriolis terms are updated.
3360 if (v_first) then
3361 is_v = isv-1 ; ie_v = iev+1 ; js_u = jsv ; je_u = jev
3362 else
3363 is_v = isv ; ie_v = iev ; js_u = jsv-1 ; je_u = jev+1
3364 endif
3365
3366 ! Use the change in eta to estimate the flow divergence and dynamic pressure.
3367 !$OMP do
3368 do j=jsv-1,jev+1 ; do i=isv-1,iev+1
3369 p_surf_dyn(i,j) = dyn_coef_eta(i,j) * (eta_pred(i,j) - eta(i,j))
3370 enddo ; enddo
3371
3372 !$OMP do schedule(static)
3373 do j=js_u,je_u ; do i=isv-1,iev
3374 pfu(i,j) = pfu(i,j) + (p_surf_dyn(i,j) - p_surf_dyn(i+1,j)) * cs%IdxCu(i,j)
3375 enddo ; enddo
3376 !$OMP end do nowait
3377 !$OMP do schedule(static)
3378 do j=jsv-1,jev ; do i=is_v,ie_v
3379 pfv(i,j) = pfv(i,j) + (p_surf_dyn(i,j) - p_surf_dyn(i,j+1)) * cs%IdyCv(i,j)
3380 enddo ; enddo
3381 !$OMP end do nowait
3382
3383end subroutine btloop_add_dyn_pf
3384
3385!> Update meridional velocity.
3386subroutine btloop_update_v(dtbt, ubt, vbt, v_accel_bt, &
3387 Cor_v, PFv, is_v, ie_v, Js_v, Je_v, f_4_v, &
3388 bt_rem_v, BT_force_v, Cor_ref_v, Rayleigh_v, &
3389 wt_accel_n, G, US, CS, Cor_bracket_bug)
3390 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3391 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3392 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3393 ubt !< The zonal barotropic velocity [L T-1 ~> m s-1].
3394 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3395 vbt !< The meridional barotropic velocity [L T-1 ~> m s-1].
3396 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3397 v_accel_bt !< The difference between the meridional acceleration from the
3398 !! barotropic calculation and BT_force_v [L T-2 ~> m s-2].
3399 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(inout) :: &
3400 Cor_v !< The meridional Coriolis acceleration [L T-2 ~> m s-2]
3401 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3402 PFv !< The meridional pressure force acceleration [L T-2 ~> m s-2].
3403 real, dimension(4,SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3404 f_4_v !< The terms giving the contribution to the Coriolis acceleration at a meridional
3405 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
3406 !! These are the products of thicknesses at u points and appropriately staggered
3407 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
3408 !! they are approximately f / 4. The 4 values on the innermost loop are for
3409 !! u-velocities to the southwest, southeast, northwest and northeast.
3410 integer, intent(in) :: is_v !< The starting i-index of the range of v-point values to calculate
3411 integer, intent(in) :: ie_v !< The ending i-index of the range of v-point values to calculate
3412 integer, intent(in) :: Js_v !< The starting j-index of the range of v-point values to calculate
3413 integer, intent(in) :: Je_v !< The ending j-index of the range of v-point values to calculate
3414 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3415 bt_rem_v !< The fraction of the barotropic meridional velocity that
3416 !! remains after a time step, the rest being lost to bottom
3417 !! drag [nondim]. bt_rem_v is between 0 and 1.
3418 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3419 BT_force_v !< The vertical average of all of the v-accelerations that are
3420 !! not explicitly included in the barotropic equation [L T-2 ~> m s-2].
3421 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3422 Cor_ref_v !< The meridional barotropic Coriolis acceleration due
3423 !! to the reference velocities [L T-2 ~> m s-2].
3424 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3425 Rayleigh_v !< A Rayleigh drag timescale operating at v-points for drag parameterizations
3426 !! that introduced directly into the barotropic solver rather than coming
3427 !! in via the visc_rem_v arrays from the layered equations [T-1 ~> s-1]
3428 real, intent(in) :: wt_accel_n !< The raw or relative weights of each of the barotropic timesteps
3429 !! in determining the average accelerations [nondim]
3430 real, intent(in) :: dtbt !< The barotropic time step [T ~> s].
3431 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3432 logical, optional, intent(in) :: Cor_bracket_bug !< If present and true, use an order of operations that is
3433 !! not bitwise rotationally symmetric in the meridional Coriolis term
3434
3435 ! Local variables
3436 logical :: use_bracket_bug
3437 integer :: i, j
3438
3439 use_bracket_bug = .false. ; if (present(cor_bracket_bug)) use_bracket_bug = cor_bracket_bug
3440
3441 ! The bracket bug only applies if v is second, use ioff to check.
3442 if (use_bracket_bug) then
3443 !$OMP do schedule(static)
3444 do j=js_v,je_v ; do i=is_v,ie_v
3445 cor_v(i,j) = -1.0*(((f_4_v(1,i,j) * ubt(i-1,j)) + (f_4_v(2,i,j) * ubt(i,j))) + &
3446 ((f_4_v(4,i,j) * ubt(i,j+1)) + (f_4_v(3,i,j) * ubt(i-1,j+1)))) - cor_ref_v(i,j)
3447 enddo ; enddo
3448 !$OMP end do nowait
3449 else
3450 !$OMP do schedule(static)
3451 do j=js_v,je_v ; do i=is_v,ie_v
3452 cor_v(i,j) = -1.0*(((f_4_v(1,i,j) * ubt(i-1,j)) + (f_4_v(4,i,j) * ubt(i,j+1))) + &
3453 ((f_4_v(2,i,j) * ubt(i,j)) + (f_4_v(3,i,j) * ubt(i-1,j+1)))) - cor_ref_v(i,j)
3454 enddo ; enddo
3455 !$OMP end do nowait
3456 endif
3457
3458 !$OMP do schedule(static)
3459 ! This updates the v-velocity, except at OBC points.
3460 do j=js_v,je_v ; do i=is_v,ie_v
3461 vbt(i,j) = bt_rem_v(i,j) * (vbt(i,j) + &
3462 dtbt * ((bt_force_v(i,j) + cor_v(i,j)) + pfv(i,j)))
3463 if (abs(vbt(i,j)) < cs%vel_underflow) vbt(i,j) = 0.0
3464 enddo ; enddo
3465 !$OMP end do nowait
3466
3467 if (cs%linear_wave_drag) then
3468 !$OMP do schedule(static)
3469 do j=js_v,je_v ; do i=is_v,ie_v
3470 v_accel_bt(i,j) = v_accel_bt(i,j) + wt_accel_n * &
3471 ((cor_v(i,j) + pfv(i,j)) - vbt(i,j)*rayleigh_v(i,j))
3472 enddo ; enddo
3473 else
3474 !$OMP do schedule(static)
3475 do j=js_v,je_v ; do i=is_v,ie_v
3476 v_accel_bt(i,j) = v_accel_bt(i,j) + wt_accel_n * (cor_v(i,j) + pfv(i,j))
3477 enddo ; enddo
3478 endif
3479
3480end subroutine btloop_update_v
3481
3482!> Update zonal velocity.
3483subroutine btloop_update_u(dtbt, ubt, vbt, u_accel_bt, &
3484 Cor_u, PFu, Is_u, Ie_u, js_u, je_u, f_4_u, &
3485 bt_rem_u, BT_force_u, Cor_ref_u, Rayleigh_u, &
3486 wt_accel_n, G, US, CS)
3487 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3488 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3489 real, intent(in) :: dtbt !< The barotropic time step [T ~> s].
3490 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3491 ubt !< The zonal barotropic velocity [L T-1 ~> m s-1].
3492 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3493 vbt !< The meridional barotropic velocity [L T-1 ~> m s-1].
3494 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3495 u_accel_bt !! The difference between the zonal acceleration from the
3496 !< barotropic calculation and BT_force_v [L T-2 ~> m s-2].
3497 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(inout) :: &
3498 Cor_u !< The anomalous zonal Coriolis acceleration [L T-2 ~> m s-2]
3499 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3500 PFu !< The anomalous zonal pressure force acceleration [L T-2 ~> m s-2].
3501 integer, intent(in) :: Is_u !< The starting i-index of the range of u-point values to calculate
3502 integer, intent(in) :: Ie_u !< The ending i-index of the range of u-point values to calculate
3503 integer, intent(in) :: js_u !< The starting j-index of the range of u-point values to calculate
3504 integer, intent(in) :: je_u !< The ending j-index of the range of u-point values to calculate
3505 real, dimension(4,SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3506 f_4_u !< The terms giving the contribution to the Coriolis acceleration at a zonal
3507 !! velocity point from the neighboring meridional velocity anomalies [T-1 ~> s-1].
3508 !! These are the products of thicknesses at v points and appropriately staggered
3509 !! averaged pseudo potential vorticities, but with sufficiently smooth topography
3510 !! they are approximately f / 4. The 4 values on the innermost loop are for
3511 !! v-velocities to the southwest, southeast, northwest and northeast.
3512 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3513 bt_rem_u !< The fraction of the barotropic meridional velocity that
3514 !! remains after a time step, the rest being lost to bottom
3515 !! drag [nondim]. bt_rem_v is between 0 and 1.
3516 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3517 BT_force_u !< The vertical average of all of the v-accelerations that are
3518 !! not explicitly included in the barotropic equation [L T-2 ~> m s-2].
3519 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3520 Cor_ref_u !< The meridional barotropic Coriolis acceleration due
3521 !! to the reference velocities [L T-2 ~> m s-2].
3522 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3523 Rayleigh_u !< A Rayleigh drag timescale operating at u-points for drag parameterizations
3524 !! that introduced directly into the barotropic solver rather than coming
3525 !! in via the visc_rem_u arrays from the layered equations [T-1 ~> s-1].
3526 real, intent(in) :: wt_accel_n !< The raw or relative weights of each of the barotropic timesteps
3527 !! in determining the average accelerations [nondim]
3528 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3529
3530 ! Local variables
3531 integer :: i, j
3532
3533 !$OMP do schedule(static)
3534 do j=js_u,je_u ; do i=is_u,ie_u
3535 cor_u(i,j) = (((f_4_u(4,i,j) * vbt(i+1,j)) + (f_4_u(1,i,j) * vbt(i,j-1))) + &
3536 ((f_4_u(3,i,j) * vbt(i,j)) + (f_4_u(2,i,j) * vbt(i+1,j-1)))) - &
3537 cor_ref_u(i,j)
3538
3539 ubt(i,j) = bt_rem_u(i,j) * (ubt(i,j) + &
3540 dtbt * ((bt_force_u(i,j) + cor_u(i,j)) + pfu(i,j)))
3541 if (abs(ubt(i,j)) < cs%vel_underflow) ubt(i,j) = 0.0
3542 enddo ; enddo
3543 !$OMP end do nowait
3544
3545 if (cs%linear_wave_drag) then
3546 !$OMP do schedule(static)
3547 do j=js_u,je_u ; do i=is_u,ie_u
3548 u_accel_bt(i,j) = u_accel_bt(i,j) + wt_accel_n * &
3549 ((cor_u(i,j) + pfu(i,j)) - ubt(i,j)*rayleigh_u(i,j))
3550 enddo ; enddo
3551 !$OMP end do nowait
3552 else
3553 !$OMP do schedule(static)
3554 do j=js_u,je_u ; do i=is_u,ie_u
3555 u_accel_bt(i,j) = u_accel_bt(i,j) + wt_accel_n * (cor_u(i,j) + pfu(i,j))
3556 enddo ; enddo
3557 !$OMP end do nowait
3558 endif
3559
3560end subroutine btloop_update_u
3561
3562
3563!> Calculate the zonal and meridional velocity from the 3-D velocity.
3564subroutine btstep_ubt_from_layer(U_in, V_in, wt_u, wt_v, ubt, vbt, G, GV, CS)
3565 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
3566 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3567 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3568 real, intent(in) :: U_in(SZIB_(G),SZJ_(G),SZK_(GV)) !< The initial (3-D) zonal velocity [L T-1 ~> m s-1]
3569 real, intent(in) :: V_in(SZI_(G),SZJB_(G),SZK_(GV)) !< The initial (3-D) meridional velocity [L T-1 ~> m s-1]
3570 real, intent(in) :: wt_u(SZIB_(G),SZJ_(G),SZK_(GV)) !< The normalized weights to be used in calculating
3571 !! zonal barotropic velocities, possibly with sums
3572 !! less than one due to viscous losses [nondim]
3573 real, intent(in) :: wt_v(SZI_(G),SZJB_(G),SZK_(GV)) !< The normalized weights to be used in calculating
3574 !! meridional barotropic velocities, possibly with
3575 !! sums less than one due to viscous losses [nondim]
3576 real, intent(out) :: ubt(SZIBW_(CS),SZJW_(CS)) !< The zonal barotropic velocity [L T-1 ~> m s-1]
3577 real, intent(out) :: vbt(SZIW_(CS),SZJBW_(CS)) !< The meridional barotropic velocity [L T-1 ~> m s-1]
3578
3579 ! Local variables
3580 integer :: i, j, k, is, ie, js, je, nz
3581
3582 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
3583
3584 ubt(:,:) = 0.0 ; vbt(:,:) = 0.0
3585
3586 !$OMP parallel do default(shared)
3587 do j=js,je ; do k=1,nz ; do i=is-1,ie
3588 ubt(i,j) = ubt(i,j) + wt_u(i,j,k) * u_in(i,j,k)
3589 enddo ; enddo ; enddo
3590 !$OMP parallel do default(shared)
3591 do j=js-1,je ; do k=1,nz ; do i=is,ie
3592 vbt(i,j) = vbt(i,j) + wt_v(i,j,k) * v_in(i,j,k)
3593 enddo ; enddo ; enddo
3594
3595 !$OMP parallel do default(shared)
3596 do j=js,je ; do i=is-1,ie
3597 if (abs(ubt(i,j)) < cs%vel_underflow) ubt(i,j) = 0.0
3598 enddo ; enddo
3599 !$OMP parallel do default(shared)
3600 do j=js-1,je ; do i=is,ie
3601 if (abs(vbt(i,j)) < cs%vel_underflow) vbt(i,j) = 0.0
3602 enddo ; enddo
3603
3604end subroutine btstep_ubt_from_layer
3605
3606
3607!> Calculate the zonal and meridional acceleration of each layer due to the barotropic calculation.
3608subroutine btstep_layer_accel(dt, u_accel_bt, v_accel_bt, pbce, gtot_E, gtot_W, gtot_N, gtot_S, &
3609 e_anom, G, GV, CS, accel_layer_u, accel_layer_v)
3610 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
3611 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
3612 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
3613 real, intent(in) :: dt !< The time increment to integrate over [T ~> s].
3614 real, dimension(SZIBW_(CS),SZJW_(CS)), intent(in) :: &
3615 u_accel_bt !< The difference between the zonal acceleration from the
3616 !! barotropic calculation and BT_force_u [L T-2 ~> m s-2].
3617 real, dimension(SZIW_(CS),SZJBW_(CS)), intent(in) :: &
3618 v_accel_bt !< The difference between the meridional acceleration from the
3619 !! barotropic calculation and BT_force_v [L T-2 ~> m s-2].
3620 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: pbce !< The baroclinic pressure anomaly in each layer
3621 !! due to free surface height anomalies
3622 !! [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3623 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3624 gtot_E !< The effective total reduced gravity used to relate free surface height
3625 !! deviations to pressure forces (including GFS and baroclinic contributions)
3626 !! in the barotropic momentum equations half a grid-point to the east of a
3627 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3628 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3629 gtot_W !< The effective total reduced gravity used to relate free surface height
3630 !! deviations to pressure forces (including GFS and baroclinic contributions)
3631 !! in the barotropic momentum equations half a grid-point to the west of a
3632 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3633 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3634 gtot_N !< The effective total reduced gravity used to relate free surface height
3635 !! deviations to pressure forces (including GFS and baroclinic contributions)
3636 !! in the barotropic momentum equations half a grid-point to the north of a
3637 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3638 real, dimension(SZIW_(CS),SZJW_(CS)), intent(in) :: &
3639 gtot_S !< The effective total reduced gravity used to relate free surface height
3640 !! deviations to pressure forces (including GFS and baroclinic contributions)
3641 !! in the barotropic momentum equations half a grid-point to the south of a
3642 !! thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3643 !! (See Hallberg, J Comp Phys 1997 for a discussion of gtot_E, etc.)
3644 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: &
3645 e_anom !< The anomaly in the sea surface height or column mass
3646 !! averaged between the beginning and end of the time step,
3647 !! relative to eta_PF, with SAL effects included [H ~> m or kg m-2].
3648 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(out) :: accel_layer_u !< The zonal acceleration of each layer due
3649 !! to the barotropic calculation [L T-2 ~> m s-2].
3650 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(out) :: accel_layer_v !< The meridional acceleration of each layer
3651 !! due to the barotropic calculation [L T-2 ~> m s-2].
3652
3653 ! Local variables
3654 real :: accel_underflow ! An acceleration that is so small it should be zeroed out [L T-2 ~> m s-2].
3655 real :: Idt ! The inverse of dt [T-1 ~> s-1].
3656 integer :: i, j, k, is, ie, js, je, nz
3657
3658 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
3659
3660 idt = 1.0 / dt
3661 accel_underflow = cs%vel_underflow * idt
3662
3663 ! Now calculate each layer's accelerations.
3664 !$OMP parallel do default(shared)
3665 do k=1,nz
3666 do j=js,je ; do i=is-1,ie
3667 accel_layer_u(i,j,k) = (u_accel_bt(i,j) - &
3668 (((pbce(i+1,j,k) - gtot_w(i+1,j)) * e_anom(i+1,j)) - &
3669 ((pbce(i,j,k) - gtot_e(i,j)) * e_anom(i,j))) * cs%IdxCu(i,j) )
3670 if (abs(accel_layer_u(i,j,k)) < accel_underflow) accel_layer_u(i,j,k) = 0.0
3671 enddo ; enddo
3672 do j=js-1,je ; do i=is,ie
3673 accel_layer_v(i,j,k) = (v_accel_bt(i,j) - &
3674 (((pbce(i,j+1,k) - gtot_s(i,j+1)) * e_anom(i,j+1)) - &
3675 ((pbce(i,j,k) - gtot_n(i,j)) * e_anom(i,j))) * cs%IdyCv(i,j) )
3676 if (abs(accel_layer_v(i,j,k)) < accel_underflow) accel_layer_v(i,j,k) = 0.0
3677 enddo ; enddo
3678 enddo
3679
3680end subroutine btstep_layer_accel
3681
3682!> This subroutine automatically determines an optimal value for dtbt based on some state of the ocean. Either pbce or
3683!! gtot_est is required to calculate gravitational acceleration. Column thickness can be estimated using BT_cont, eta,
3684!! and SSH_add (default=0), with priority given in that order. The subroutine sets CS%dtbt_max and CS%dtbt.
3685subroutine set_dtbt(G, GV, US, CS, pbce, gtot_est, BT_cont, eta, SSH_add)
3686 type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
3687 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
3688 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
3689 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
3690 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
3691 optional, intent(in) :: pbce !< The baroclinic pressure anomaly in each layer due to free
3692 !! surface height anomalies [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3693 real, optional, intent(in) :: gtot_est !< An estimate of the total gravitational acceleration
3694 !! [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3695 type(bt_cont_type), optional, pointer :: bt_cont !< A structure with elements that describe the effective open
3696 !! face areas as a function of barotropic flow.
3697 real, dimension(SZI_(G),SZJ_(G)), &
3698 optional, intent(in) :: eta !< The barotropic free surface height anomaly or column mass
3699 !! anomaly [H ~> m or kg m-2].
3700 real, optional, intent(in) :: ssh_add !< An additional contribution to SSH to provide a margin of
3701 !! error when calculating the external wave speed [Z ~> m].
3702
3703 ! Local variables
3704 real, dimension(SZI_(G),SZJ_(G)) :: &
3705 gtot_e, & ! gtot_X is the effective total reduced gravity used to relate
3706 gtot_w, & ! free surface height deviations to pressure forces (including
3707 gtot_n, & ! GFS and baroclinic contributions) in the barotropic momentum
3708 gtot_s ! equations half a grid-point in the X-direction (X is N, S, E, or W)
3709 ! from the thickness point [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
3710 ! (See Hallberg, J Comp Phys 1997 for a discussion.)
3711 real, dimension(SZIBS_(G),SZJ_(G)) :: &
3712 datu ! Basin depth at u-velocity grid points times the y-grid
3713 ! spacing [H L ~> m2 or kg m-1].
3714 real, dimension(SZI_(G),SZJBS_(G)) :: &
3715 datv ! Basin depth at v-velocity grid points times the x-grid
3716 ! spacing [H L ~> m2 or kg m-1].
3717 real :: det_de ! The partial derivative due to self-attraction and loading
3718 ! of the reference geopotential with the sea surface height [nondim].
3719 ! This is typically ~0.09 or less.
3720 real :: dgeo_de ! The constant of proportionality between geopotential and
3721 ! sea surface height [nondim]. It is a nondimensional number of
3722 ! order 1. For stability, this may be made larger
3723 ! than physical problem would suggest.
3724 real :: add_ssh ! An additional contribution to SSH to provide a margin of error
3725 ! when calculating the external wave speed [Z ~> m].
3726 real :: min_max_dt2 ! The square of the minimum value of the largest stable barotropic
3727 ! timesteps [T2 ~> s2]
3728 real :: dtbt_max ! The maximum barotropic timestep [T ~> s]
3729 real :: idt_max2 ! The squared inverse of the local maximum stable
3730 ! barotropic time step [T-2 ~> s-2].
3731 logical :: use_bt_cont
3732 type(memory_size_type) :: ms
3733
3734 integer :: i, j, k, is, ie, js, je, nz
3735
3736 if (.not.cs%module_is_initialized) call mom_error(fatal, &
3737 "set_dtbt: Module MOM_barotropic must be initialized before it is used.")
3738
3739 if (.not.cs%split) return
3740 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
3741 ms%isdw = g%isd ; ms%iedw = g%ied ; ms%jsdw = g%jsd ; ms%jedw = g%jed
3742
3743
3744 if (.not.(present(pbce) .or. present(gtot_est))) call mom_error(fatal, &
3745 "set_dtbt: Either pbce or gtot_est must be present.")
3746
3747 add_ssh = 0.0 ; if (present(ssh_add)) add_ssh = ssh_add
3748
3749 use_bt_cont = .false.
3750 if (present(bt_cont)) use_bt_cont = (associated(bt_cont))
3751
3752 if (use_bt_cont) then
3753 call bt_cont_to_face_areas(bt_cont, datu, datv, g, us, ms, halo=0)
3754 elseif (cs%Nonlinear_continuity .and. present(eta)) then
3755 call find_face_areas(datu, datv, g, gv, us, cs, ms, 0, eta=eta)
3756 else
3757 call find_face_areas(datu, datv, g, gv, us, cs, ms, 0, add_max=add_ssh)
3758 endif
3759
3760 det_de = 0.0
3761 if (cs%calculate_SAL) call scalar_sal_sensitivity(cs%SAL_CSp, det_de)
3762 if (cs%tidal_sal_bug) then
3763 dgeo_de = 1.0 + max(0.0, det_de + cs%G_extra)
3764 else
3765 dgeo_de = 1.0 + max(0.0, cs%G_extra - det_de)
3766 endif
3767 if (present(pbce)) then
3768 do j=js,je ; do i=is,ie
3769 gtot_e(i,j) = 0.0 ; gtot_w(i,j) = 0.0
3770 gtot_n(i,j) = 0.0 ; gtot_s(i,j) = 0.0
3771 enddo ; enddo
3772 do k=1,nz ; do j=js,je ; do i=is,ie
3773 gtot_e(i,j) = gtot_e(i,j) + pbce(i,j,k) * cs%frhatu(i,j,k)
3774 gtot_w(i,j) = gtot_w(i,j) + pbce(i,j,k) * cs%frhatu(i-1,j,k)
3775 gtot_n(i,j) = gtot_n(i,j) + pbce(i,j,k) * cs%frhatv(i,j,k)
3776 gtot_s(i,j) = gtot_s(i,j) + pbce(i,j,k) * cs%frhatv(i,j-1,k)
3777 enddo ; enddo ; enddo
3778 else
3779 do j=js,je ; do i=is,ie
3780 gtot_e(i,j) = gtot_est ; gtot_w(i,j) = gtot_est
3781 gtot_n(i,j) = gtot_est ; gtot_s(i,j) = gtot_est
3782 enddo ; enddo
3783 endif
3784
3785 min_max_dt2 = 1.0e38*us%s_to_T**2 ! A huge value for the permissible timestep squared.
3786 do j=js,je ; do i=is,ie
3787 ! This is pretty accurate for gravity waves, but it is a conservative
3788 ! estimate since it ignores the stabilizing effect of the bottom drag.
3789 idt_max2 = 0.5 * (1.0 + 2.0*cs%bebt) * (g%IareaT(i,j) * &
3790 (((gtot_e(i,j)*datu(i,j)*g%IdxCu(i,j)) + (gtot_w(i,j)*datu(i-1,j)*g%IdxCu(i-1,j))) + &
3791 ((gtot_n(i,j)*datv(i,j)*g%IdyCv(i,j)) + (gtot_s(i,j)*datv(i,j-1)*g%IdyCv(i,j-1)))) + &
3792 ((g%Coriolis2Bu(i,j) + g%Coriolis2Bu(i-1,j-1)) + &
3793 (g%Coriolis2Bu(i-1,j) + g%Coriolis2Bu(i,j-1))) * cs%BT_Coriolis_scale**2 )
3794 if (idt_max2 * min_max_dt2 > 1.0) min_max_dt2 = 1.0 / idt_max2
3795 enddo ; enddo
3796 dtbt_max = sqrt(min_max_dt2 / dgeo_de)
3797 if (id_clock_sync > 0) call cpu_clock_begin(id_clock_sync)
3798 call min_across_pes(dtbt_max)
3799 if (id_clock_sync > 0) call cpu_clock_end(id_clock_sync)
3800
3801 cs%dtbt = cs%dtbt_fraction * dtbt_max
3802 cs%dtbt_max = dtbt_max
3803
3804 if (cs%debug) then
3805 call chksum0(cs%dtbt, "End set_dtbt dtbt", unscale=us%T_to_s)
3806 call chksum0(cs%dtbt_max, "End set_dtbt dtbt_max", unscale=us%T_to_s)
3807 endif
3808
3809end subroutine set_dtbt
3810
3811! The following 5 subroutines apply the open boundary conditions.
3812
3813!> This subroutine applies the open boundary conditions on barotropic zonal
3814!! velocities and mass transports, as developed by Mehmet Ilicak.
3815subroutine apply_u_velocity_obcs(ubt, uhbt, ubt_trans, eta, SpV_avg, ubt_old, BT_OBC, G, MS, &
3816 GV, US, CS, halo, dtbt, bebt, use_BT_cont, integral_BT_cont, dt_elapsed, &
3817 Datu, BTCL_u, uhbt0, ubt_int, ubt_int_prev, uhbt_int, uhbt_int_prev)
3818 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
3819 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of
3820 !! the argument arrays.
3821 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(inout) :: ubt !< the zonal barotropic velocity [L T-1 ~> m s-1].
3822 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(inout) :: uhbt !< the zonal barotropic transport
3823 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
3824 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(inout) :: ubt_trans !< The zonal barotropic velocity used in
3825 !! transport [L T-1 ~> m s-1].
3826 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: eta !< The barotropic free surface height anomaly or
3827 !! column mass anomaly [H ~> m or kg m-2].
3828 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: SpV_avg !< The column average specific volume [R-1 ~> m3 kg-1]
3829 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: ubt_old !< The starting value of ubt in a barotropic
3830 !! step [L T-1 ~> m s-1].
3831 type(bt_obc_type), intent(in) :: BT_OBC !< A structure with the private barotropic arrays
3832 !! related to the open boundary conditions,
3833 !! set by set_up_BT_OBC.
3834 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
3835 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
3836 type(barotropic_cs), intent(in) :: CS !< Barotropic control structure
3837 integer, intent(in) :: halo !< The extra halo size to use here.
3838 real, intent(in) :: dtbt !< The time step [T ~> s].
3839 real, intent(in) :: bebt !< The fractional weighting of the future velocity
3840 !! in determining the transport [nondim]
3841 logical, intent(in) :: use_BT_cont !< If true, use the BT_cont_types to calculate
3842 !! transports.
3843 logical, intent(in) :: integral_BT_cont !< If true, update the barotropic continuity
3844 !! equation directly from the initial condition
3845 !! using the time-integrated barotropic velocity.
3846 real, intent(in) :: dt_elapsed !< The amount of time in the barotropic stepping
3847 !! that will have elapsed [T ~> s].
3848 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: Datu !< A fixed estimate of the face areas at u points
3849 !! [H L ~> m2 or kg m-1].
3850 type(local_bt_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: BTCL_u !< Structure of information used
3851 !! for a dynamic estimate of the face areas at
3852 !! u-points.
3853 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: uhbt0 !< A correction to the zonal transport so that
3854 !! the barotropic functions agree with the sum
3855 !! of the layer transports
3856 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
3857 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(inout) :: ubt_int !< The time-integrated zonal barotropic
3858 !! velocity after this update [L T-1 ~> m s-1]
3859 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: ubt_int_prev !< The time-integrated zonal barotropic
3860 !! velocity before this update [L T-1 ~> m s-1]
3861 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(inout) :: uhbt_int !< The time-integrated zonal barotropic transport
3862 !! after this update [H L2 T-1 ~> m3 s-1 or kg s-1]
3863 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: uhbt_int_prev !< The time-integrated zonal barotropic
3864 !! transport before this update
3865 !! [H L2 T-1 ~> m3 s-1 or kg s-1]
3866
3867 ! Local variables
3868 real :: vel_prev ! The previous velocity [L T-1 ~> m s-1].
3869 real :: cfl ! The CFL number at the point in question [nondim]
3870 real :: u_inlet ! The zonal inflow velocity [L T-1 ~> m s-1]
3871 real :: uhbt_int_new ! The updated time-integrated zonal transport [H L2 ~> m3]
3872 real :: ssh_in ! The inflow sea surface height [Z ~> m]
3873 real :: ssh_1 ! The sea surface height in the interior cell adjacent to the an OBC face [Z ~> m]
3874 real :: ssh_2 ! The sea surface height in the next cell inward from the OBC face [Z ~> m]
3875 real :: Idtbt ! The inverse of the barotropic time step [T-1 ~> s-1]
3876 integer :: i, j, Is_u, Ie_u, js, je
3877
3878 if (.not.bt_obc%u_OBCs_on_PE) return
3879
3880 idtbt = 1.0 / dtbt
3881
3882 ! Work on Eastern OBC points
3883 is_u = max((g%isc-1)-halo, bt_obc%Is_u_E_obc) ; ie_u = min(g%iec+halo, bt_obc%Ie_u_E_obc)
3884 js = max(g%jsc-halo, bt_obc%js_u_E_obc) ; je = min(g%jec+halo, bt_obc%je_u_E_obc)
3885 do j=js,je ; do i=is_u,ie_u ; if (bt_obc%u_OBC_type(i,j) > 0) then
3886 if (bt_obc%u_OBC_type(i,j) == specified_obc) then ! Eastern specified OBC
3887 uhbt(i,j) = bt_obc%uhbt(i,j)
3888 ubt(i,j) = bt_obc%ubt_outer(i,j)
3889 ubt_trans(i,j) = ubt(i,j)
3890 if (integral_bt_cont) then
3891 uhbt_int(i,j) = uhbt_int_prev(i,j) + dtbt * uhbt(i,j)
3892 ubt_int(i,j) = ubt_int_prev(i,j) + dtbt * ubt_trans(i,j)
3893 endif
3894 elseif (bt_obc%u_OBC_type(i,j) == flather_obc) then ! Eastern Flather OBC
3895 cfl = dtbt * bt_obc%Cg_u(i,j) * g%IdxCu(i,j) ! CFL
3896 u_inlet = cfl*ubt_old(i-1,j) + (1.0-cfl)*ubt_old(i,j) ! Valid for cfl<1
3897 if (i <= ms%isdw) then
3898 ! Do not apply an Eastern Flather OBC at the western halo points on a PE, as doing so would
3899 ! create a segmentation fault and this velocity will not propagate through to the next iteration.
3900 ssh_in = bt_obc%SSH_outer_u(i,j)
3901 elseif (gv%Boussinesq) then
3902 ssh_in = gv%H_to_Z*(eta(i,j) + (0.5-cfl)*(eta(i,j)-eta(i-1,j))) ! internal
3903 else
3904 ssh_1 = gv%H_to_RZ * eta(i,j) * spv_avg(i,j) - (cs%bathyT(i,j) + g%Z_ref)
3905 ssh_2 = gv%H_to_RZ * eta(i-1,j) * spv_avg(i-1,j) - (cs%bathyT(i-1,j) + g%Z_ref)
3906 ssh_in = ssh_1 + (0.5-cfl)*(ssh_1-ssh_2) ! internal
3907 endif
3908 if (bt_obc%dZ_u(i,j) > 0.0) then
3909 vel_prev = ubt(i,j)
3910 ubt(i,j) = 0.5*((u_inlet + bt_obc%ubt_outer(i,j)) + &
3911 (bt_obc%Cg_u(i,j)/bt_obc%dZ_u(i,j)) * (ssh_in-bt_obc%SSH_outer_u(i,j)))
3912 ubt_trans(i,j) = (1.0-bebt)*vel_prev + bebt*ubt(i,j)
3913 else ! This point is now dry.
3914 ubt(i,j) = 0.0
3915 ubt_trans(i,j) = 0.0
3916 endif
3917 elseif (bt_obc%u_OBC_type(i,j) == gradient_obc) then ! Eastern gradient OBC
3918 ubt(i,j) = ubt(i-1,j)
3919 ubt_trans(i,j) = ubt(i,j)
3920 endif
3921
3922 ! Reset transports and related time-inetegrated velocities with non-specified OBCs
3923 if (bt_obc%u_OBC_type(i,j) > specified_obc) then ! Eastern Flather or gradient OBC
3924 if (integral_bt_cont) then
3925 ubt_int(i,j) = ubt_int_prev(i,j) + dtbt * ubt_trans(i,j)
3926 uhbt_int_new = find_uhbt(ubt_int(i,j), btcl_u(i,j)) + dt_elapsed*uhbt0(i,j)
3927 uhbt(i,j) = (uhbt_int_new - uhbt_int_prev(i,j)) * idtbt
3928 uhbt_int(i,j) = uhbt_int_prev(i,j) + dtbt * uhbt(i,j)
3929 ! The line above is equivalent to: uhbt_int(I,j) = uhbt_int_new
3930 elseif (use_bt_cont) then
3931 uhbt(i,j) = find_uhbt(ubt_trans(i,j), btcl_u(i,j)) + uhbt0(i,j)
3932 else
3933 uhbt(i,j) = datu(i,j)*ubt_trans(i,j) + uhbt0(i,j)
3934 endif
3935 endif
3936
3937 endif ; enddo ; enddo
3938
3939 ! Work on Western OBC points
3940 is_u = max((g%isc-1)-halo, bt_obc%Is_u_W_obc) ; ie_u = min(g%iec+halo, bt_obc%Ie_u_W_obc)
3941 js = max(g%jsc-halo, bt_obc%js_u_W_obc) ; je = min(g%jec+halo, bt_obc%je_u_W_obc)
3942 do j=js,je ; do i=is_u,ie_u ; if (bt_obc%u_OBC_type(i,j) < 0) then
3943 if (bt_obc%u_OBC_type(i,j) == -specified_obc) then ! Western specified OBC
3944 uhbt(i,j) = bt_obc%uhbt(i,j)
3945 ubt(i,j) = bt_obc%ubt_outer(i,j)
3946 ubt_trans(i,j) = ubt(i,j)
3947 if (integral_bt_cont) then
3948 uhbt_int(i,j) = uhbt_int_prev(i,j) + dtbt * uhbt(i,j)
3949 ubt_int(i,j) = ubt_int_prev(i,j) + dtbt * ubt_trans(i,j)
3950 endif
3951 elseif (bt_obc%u_OBC_type(i,j) == -flather_obc) then ! Western Flather OBC
3952 cfl = dtbt * bt_obc%Cg_u(i,j) * g%IdxCu(i,j) ! CFL
3953 u_inlet = cfl*ubt_old(i+1,j) + (1.0-cfl)*ubt_old(i,j) ! Valid for cfl<1
3954 if (i >= ms%iedw-1) then
3955 ! Do not apply a Western Flather OBC at the eastern halo points on a PE, as doing so would
3956 ! create a segmentation fault and this velocity will not propagate through to the next iteration.
3957 ssh_in = bt_obc%SSH_outer_u(i,j)
3958 elseif (gv%Boussinesq) then
3959 ssh_in = gv%H_to_Z*(eta(i+1,j) + (0.5-cfl)*(eta(i+1,j)-eta(i+2,j))) ! internal
3960 else
3961 ssh_1 = gv%H_to_RZ * eta(i+1,j) * spv_avg(i+1,j) - (cs%bathyT(i+1,j) + g%Z_ref)
3962 ssh_2 = gv%H_to_RZ * eta(i+2,j) * spv_avg(i+2,j) - (cs%bathyT(i+2,j) + g%Z_ref)
3963 ssh_in = ssh_1 + (0.5-cfl)*(ssh_1-ssh_2) ! internal
3964 endif
3965
3966 if (bt_obc%dZ_u(i,j) > 0.0) then
3967 vel_prev = ubt(i,j)
3968 ubt(i,j) = 0.5*((u_inlet + bt_obc%ubt_outer(i,j)) + &
3969 (bt_obc%Cg_u(i,j)/bt_obc%dZ_u(i,j)) * (bt_obc%SSH_outer_u(i,j)-ssh_in))
3970 ubt_trans(i,j) = (1.0-bebt)*vel_prev + bebt*ubt(i,j)
3971 else ! This point is now dry.
3972 ubt(i,j) = 0.0
3973 ubt_trans(i,j) = 0.0
3974 endif
3975 elseif (bt_obc%u_OBC_type(i,j) == -gradient_obc) then ! Western gradient OBC
3976 ubt(i,j) = ubt(i+1,j)
3977 ubt_trans(i,j) = ubt(i,j)
3978 endif
3979
3980 ! Reset transports and related time-inetegrated velocities with non-specified OBCs
3981 if (bt_obc%u_OBC_type(i,j) < -specified_obc) then ! Western Flather or gradient OBC
3982 if (integral_bt_cont) then
3983 ubt_int(i,j) = ubt_int_prev(i,j) + dtbt * ubt_trans(i,j)
3984 uhbt_int_new = find_uhbt(ubt_int(i,j), btcl_u(i,j)) + dt_elapsed*uhbt0(i,j)
3985 uhbt(i,j) = (uhbt_int_new - uhbt_int_prev(i,j)) * idtbt
3986 uhbt_int(i,j) = uhbt_int_prev(i,j) + dtbt * uhbt(i,j)
3987 ! The line above is equivalent to: uhbt_int(I,j) = uhbt_int_new
3988 elseif (use_bt_cont) then
3989 uhbt(i,j) = find_uhbt(ubt_trans(i,j), btcl_u(i,j)) + uhbt0(i,j)
3990 else
3991 uhbt(i,j) = datu(i,j)*ubt_trans(i,j) + uhbt0(i,j)
3992 endif
3993 endif
3994
3995 endif ; enddo ; enddo
3996
3997end subroutine apply_u_velocity_obcs
3998
3999!> This subroutine applies the open boundary conditions on barotropic meridional
4000!! velocities and mass transports, as developed by Mehmet Ilicak.
4001subroutine apply_v_velocity_obcs(vbt, vhbt, vbt_trans, eta, SpV_avg, vbt_old, BT_OBC, &
4002 G, MS, GV, US, CS, halo, dtbt, bebt, use_BT_cont, integral_BT_cont, dt_elapsed, &
4003 Datv, BTCL_v, vhbt0, vbt_int, vbt_int_prev, vhbt_int, vhbt_int_prev)
4004 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
4005 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of
4006 !! the argument arrays.
4007 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(inout) :: vbt !< The meridional barotropic velocity
4008 !! [L T-1 ~> m s-1].
4009 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(inout) :: vhbt !< the meridional barotropic transport
4010 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
4011 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(inout) :: vbt_trans !< the meridional BT velocity used in
4012 !! transports [L T-1 ~> m s-1].
4013 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: eta !< The barotropic free surface height anomaly or
4014 !! column mass anomaly [H ~> m or kg m-2].
4015 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: SpV_avg !< The column average specific volume [R-1 ~> m3 kg-1]
4016 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vbt_old !< The starting value of vbt in a barotropic
4017 !! step [L T-1 ~> m s-1].
4018 type(bt_obc_type), intent(in) :: BT_OBC !< A structure with the private barotropic arrays
4019 !! related to the open boundary conditions,
4020 !! set by set_up_BT_OBC.
4021 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
4022 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
4023 type(barotropic_cs), intent(in) :: CS !< Barotropic control structure
4024 integer, intent(in) :: halo !< The extra halo size to use here.
4025 real, intent(in) :: dtbt !< The time step [T ~> s].
4026 real, intent(in) :: bebt !< The fractional weighting of the future velocity
4027 !! in determining the transport [nondim]
4028 logical, intent(in) :: use_BT_cont !< If true, use the BT_cont_types to calculate
4029 !! transports.
4030 logical, intent(in) :: integral_BT_cont !< If true, update the barotropic continuity
4031 !! equation directly from the initial condition
4032 !! using the time-integrated barotropic velocity.
4033 real, intent(in) :: dt_elapsed !< The amount of time in the barotropic stepping
4034 !! that will have elapsed [T ~> s].
4035 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: Datv !< A fixed estimate of the face areas at v points
4036 !! [H L ~> m2 or kg m-1].
4037 type(local_bt_cont_v_type), dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: BTCL_v !< Structure of information used
4038 !! for a dynamic estimate of the face areas at
4039 !! v-points.
4040 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vhbt0 !< A correction to the meridional transport so that
4041 !! the barotropic functions agree with the sum
4042 !! of the layer transports
4043 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
4044 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(inout) :: vbt_int !< The time-integrated meridional barotropic
4045 !! velocity after this update [L T-1 ~> m s-1].
4046 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vbt_int_prev !< The time-integrated meridional barotropic
4047 !! velocity before this update [L T-1 ~> m s-1].
4048 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(inout) :: vhbt_int !< The time-integrated meridional barotropic
4049 !! transport after this update
4050 !! [H L2 T-1 ~> m3 s-1 or kg s-1]
4051 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vhbt_int_prev !< The time-integrated meridional barotropic
4052 !! transport before this update
4053 !! [H L2 T-1 ~> m3 s-1 or kg s-1]
4054
4055 ! Local variables
4056 real :: vel_prev ! The previous velocity [L T-1 ~> m s-1].
4057 real :: cfl ! The CFL number at the point in question [nondim]
4058 real :: v_inlet ! The meridional inflow velocity [L T-1 ~> m s-1]
4059 real :: vhbt_int_new ! The updated time-integrated meridional transport [H L2 ~> m3]
4060 real :: ssh_in ! The inflow sea surface height [Z ~> m]
4061 real :: ssh_1 ! The sea surface height in the interior cell adjacent to the an OBC face [Z ~> m]
4062 real :: ssh_2 ! The sea surface height in the next cell inward from the OBC face [Z ~> m]
4063 real :: Idtbt ! The inverse of the barotropic time step [T-1 ~> s-1]
4064 integer :: i, j, is, ie, Js_v, Je_v
4065
4066 if (.not.bt_obc%v_OBCs_on_PE) return
4067
4068 idtbt = 1.0 / dtbt
4069
4070 ! This routine uses separate blocks of code and loops for Northern and southern open boundary
4071 ! condition points, despite this leading to some code duplication, because the OBCs almost always
4072 ! occur at the edge of the domain, and in parallel appliations, most PEs will only have one or
4073 ! the other.
4074
4075
4076 ! Work on Northern OBC points
4077 is = max(g%isc-halo, bt_obc%is_v_N_obc) ; ie = min(g%iec+halo, bt_obc%ie_v_N_obc)
4078 js_v = max((g%jsc-1)-halo, bt_obc%Js_v_N_obc) ; je_v = min(g%jec+halo, bt_obc%Je_v_N_obc)
4079 do j=js_v,je_v ; do i=is,ie ; if (bt_obc%v_OBC_type(i,j) > 0) then
4080 if (bt_obc%v_OBC_type(i,j) == specified_obc) then ! Northern specified OBC
4081 vhbt(i,j) = bt_obc%vhbt(i,j)
4082 vbt(i,j) = bt_obc%vbt_outer(i,j)
4083 vbt_trans(i,j) = vbt(i,j)
4084 if (integral_bt_cont) then
4085 vbt_int(i,j) = vbt_int_prev(i,j) + dtbt * vbt(i,j)
4086 vhbt_int(i,j) = vhbt_int_prev(i,j) + dtbt * vhbt(i,j)
4087 endif
4088 elseif (bt_obc%v_OBC_type(i,j) == flather_obc) then ! Northern Flather OBC
4089 cfl = dtbt * bt_obc%Cg_v(i,j) * g%IdyCv(i,j) ! CFL
4090 v_inlet = cfl*vbt_old(i,j-1) + (1.0-cfl)*vbt_old(i,j) ! Valid for cfl<1
4091 if (j <= ms%jsdw) then
4092 ! Do not apply a Northern Flather OBC at the southern halo points on a PE, as doing so would
4093 ! create a segmentation fault and this velocity will not propagate through to the next iteration.
4094 ssh_in = bt_obc%SSH_outer_v(i,j)
4095 elseif (gv%Boussinesq) then
4096 ssh_in = gv%H_to_Z*(eta(i,j) + (0.5-cfl)*(eta(i,j)-eta(i,j-1))) ! internal
4097 else
4098 ssh_1 = gv%H_to_RZ * eta(i,j) * spv_avg(i,j) - (cs%bathyT(i,j) + g%Z_ref)
4099 ssh_2 = gv%H_to_RZ * eta(i,j-1) * spv_avg(i,j-1) - (cs%bathyT(i,j-1) + g%Z_ref)
4100 ssh_in = ssh_1 + (0.5-cfl)*(ssh_1-ssh_2) ! internal
4101 endif
4102
4103 if (bt_obc%dZ_v(i,j) > 0.0) then
4104 vel_prev = vbt(i,j)
4105 vbt(i,j) = 0.5*((v_inlet + bt_obc%vbt_outer(i,j)) + &
4106 (bt_obc%Cg_v(i,j)/bt_obc%dZ_v(i,j)) * (ssh_in-bt_obc%SSH_outer_v(i,j)))
4107 vbt_trans(i,j) = (1.0-bebt)*vel_prev + bebt*vbt(i,j)
4108 else ! This point is now dry
4109 vbt(i,j) = 0.0
4110 vbt_trans(i,j) = 0.0
4111 endif
4112 elseif (bt_obc%v_OBC_type(i,j) == gradient_obc) then ! Northern gradient OBC
4113 vbt(i,j) = vbt(i,j-1)
4114 vbt_trans(i,j) = vbt(i,j)
4115 endif
4116
4117 ! Reset transports and related time-inetegrated velocities with non-specified OBCs
4118 if (bt_obc%v_OBC_type(i,j) > specified_obc) then ! Northern Flather or gradient OBC
4119 if (integral_bt_cont) then
4120 vbt_int(i,j) = vbt_int_prev(i,j) + dtbt * vbt_trans(i,j)
4121 vhbt_int_new = find_vhbt(vbt_int(i,j), btcl_v(i,j)) + dt_elapsed*vhbt0(i,j)
4122 vhbt(i,j) = (vhbt_int_new - vhbt_int_prev(i,j)) * idtbt
4123 vhbt_int(i,j) = vhbt_int_prev(i,j) + dtbt * vhbt(i,j)
4124 ! The line above is equivalent to: vhbt_int(i,J) = vhbt_int_new
4125 elseif (use_bt_cont) then
4126 vhbt(i,j) = find_vhbt(vbt_trans(i,j), btcl_v(i,j)) + vhbt0(i,j)
4127 else
4128 vhbt(i,j) = vbt_trans(i,j)*datv(i,j) + vhbt0(i,j)
4129 endif
4130 endif
4131
4132 endif ; enddo ; enddo
4133
4134 ! Work on Southern OBC points
4135 is = max(g%isc-halo, bt_obc%is_v_S_obc) ; ie = min(g%iec+halo, bt_obc%ie_v_S_obc)
4136 js_v = max((g%jsc-1)-halo, bt_obc%Js_v_S_obc) ; je_v = min(g%jec+halo, bt_obc%Je_v_S_obc)
4137 do j=js_v,je_v ; do i=is,ie ; if (bt_obc%v_OBC_type(i,j) < 0) then
4138 if (bt_obc%v_OBC_type(i,j) == -specified_obc) then ! Southern specified OBC
4139 vhbt(i,j) = bt_obc%vhbt(i,j)
4140 vbt(i,j) = bt_obc%vbt_outer(i,j)
4141 vbt_trans(i,j) = vbt(i,j)
4142 if (integral_bt_cont) then
4143 vbt_int(i,j) = vbt_int_prev(i,j) + dtbt * vbt(i,j)
4144 vhbt_int(i,j) = vhbt_int_prev(i,j) + dtbt * vhbt(i,j)
4145 endif
4146 elseif (bt_obc%v_OBC_type(i,j) == -flather_obc) then ! Southern Flather OBC
4147 cfl = dtbt * bt_obc%Cg_v(i,j) * g%IdyCv(i,j) ! CFL
4148 v_inlet = cfl*vbt_old(i,j+1) + (1.0-cfl)*vbt_old(i,j) ! Valid for cfl <1
4149 if (j >= ms%jedw-1) then
4150 ! Do not apply a Southern Flather OBC at the northern halo points on a PE, as doing so would
4151 ! create a segmentation fault and this velocity will not propagate through to the next iteration.
4152 ssh_in = bt_obc%SSH_outer_v(i,j)
4153 elseif (gv%Boussinesq) then
4154 ssh_in = gv%H_to_Z*(eta(i,j+1) + (0.5-cfl)*(eta(i,j+1)-eta(i,j+2))) ! internal
4155 else
4156 ssh_1 = gv%H_to_RZ * eta(i,j+1) * spv_avg(i,j+1) - (cs%bathyT(i,j+1) + g%Z_ref)
4157 ssh_2 = gv%H_to_RZ * eta(i,j+2) * spv_avg(i,j+2) - (cs%bathyT(i,j+2) + g%Z_ref)
4158 ssh_in = ssh_1 + (0.5-cfl)*(ssh_1-ssh_2) ! internal
4159 endif
4160
4161 if (bt_obc%dZ_v(i,j) > 0.0) then
4162 vel_prev = vbt(i,j)
4163 vbt(i,j) = 0.5*((v_inlet + bt_obc%vbt_outer(i,j)) + &
4164 (bt_obc%Cg_v(i,j)/bt_obc%dZ_v(i,j)) * (bt_obc%SSH_outer_v(i,j)-ssh_in))
4165 vbt_trans(i,j) = (1.0-bebt)*vel_prev + bebt*vbt(i,j)
4166 else ! This point is now dry
4167 vbt(i,j) = 0.0
4168 vbt_trans(i,j) = 0.0
4169 endif
4170 elseif (bt_obc%v_OBC_type(i,j) == -gradient_obc) then ! Southern gradient OBC
4171 vbt(i,j) = vbt(i,j+1)
4172 vbt_trans(i,j) = vbt(i,j)
4173 endif
4174
4175 ! Reset transports and related time-inetegrated velocities with non-specified OBCs
4176 if (bt_obc%v_OBC_type(i,j) < -specified_obc) then ! Southern Flather or gradient OBC
4177 if (integral_bt_cont) then
4178 vbt_int(i,j) = vbt_int_prev(i,j) + dtbt * vbt_trans(i,j)
4179 vhbt_int_new = find_vhbt(vbt_int(i,j), btcl_v(i,j)) + dt_elapsed*vhbt0(i,j)
4180 vhbt(i,j) = (vhbt_int_new - vhbt_int_prev(i,j)) * idtbt
4181 vhbt_int(i,j) = vhbt_int_prev(i,j) + dtbt * vhbt(i,j)
4182 ! The line above is equivalent to: vhbt_int(i,J) = vhbt_int_new
4183 elseif (use_bt_cont) then
4184 vhbt(i,j) = find_vhbt(vbt_trans(i,j), btcl_v(i,j)) + vhbt0(i,j)
4185 else
4186 vhbt(i,j) = vbt_trans(i,j)*datv(i,j) + vhbt0(i,j)
4187 endif
4188 endif
4189
4190 endif ; enddo ; enddo
4191
4192end subroutine apply_v_velocity_obcs
4193
4194!> This subroutine sets up the time-invariant control information about the open boundary
4195!! conditions on the full wide halo domain used by the barotropic solver.
4196subroutine initialize_bt_obc(OBC, BT_OBC, G, CS)
4197 type(ocean_obc_type), target, intent(inout) :: OBC !< An associated pointer to an OBC type.
4198 type(bt_obc_type), intent(inout) :: BT_OBC !< A structure with the private barotropic arrays
4199 !! related to the open boundary conditions,
4200 !! set by set_up_BT_OBC.
4201 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
4202 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
4203
4204 ! Local variables
4205 real, dimension(SZIBW_(CS),SZJW_(CS)) :: &
4206 u_OBC ! A set of integers encoding the nature of the u-point open boundary conditions,
4207 ! converted to real numbers to work with the MOM6 halo update code [nondim]
4208 real, dimension(SZIW_(CS),SZJBW_(CS)) :: &
4209 v_OBC ! A set of integers encoding the nature of the v-point open boundary conditions,
4210 ! converted to real numbers to work with the MOM6 halo update code [nondim]
4211 integer :: OBC_type ! The integer encoding the type of OBC being used at a point [nondim]
4212 logical :: reversed_OBCs ! True of there any OBCs in the opposite halo on this PE, e.g. points
4213 ! with a southern OBC in a northern halo.
4214 logical :: any_reversed_OBCs
4215 integer :: i, j, isdw, iedw, jsdw, jedw
4216 integer :: l_seg, Flather_OBC_in_halo
4217
4218 isdw = cs%isdw ; iedw = cs%iedw ; jsdw = cs%jsdw ; jedw = cs%jedw
4219
4220 u_obc(:,:) = 0.0
4221 v_obc(:,:) = 0.0
4222
4223 do j=g%jsc,g%jec ; do i=g%isc-1,g%iec
4224
4225 obc_type = 0
4226 if (obc%segnum_u(i,j) /= 0) then
4227 l_seg = abs(obc%segnum_u(i,j))
4228 if (obc%segment(l_seg)%gradient) obc_type = gradient_obc
4229 if (obc%segment(l_seg)%Flather) obc_type = flather_obc
4230 if (obc%segment(l_seg)%specified) obc_type = specified_obc
4231 u_obc(i,j) = sign(obc_type, obc%segnum_u(i,j))
4232 endif
4233 enddo ; enddo
4234
4235 do j=g%jsc-1,g%jec ; do i=g%isc,g%iec
4236 obc_type = 0
4237 if (obc%segnum_v(i,j) /= 0) then
4238 l_seg = abs(obc%segnum_v(i,j))
4239 if (obc%segment(l_seg)%gradient) obc_type = gradient_obc
4240 if (obc%segment(l_seg)%Flather) obc_type = flather_obc
4241 if (obc%segment(l_seg)%specified) obc_type = specified_obc
4242 v_obc(i,j) = sign(obc_type, obc%segnum_v(i,j))
4243 endif
4244 enddo ; enddo
4245
4246 call pass_vector(u_obc, v_obc, cs%BT_Domain)
4247
4248 allocate(bt_obc%u_OBC_type(isdw-1:iedw,jsdw:jedw), source=0)
4249 allocate(bt_obc%v_OBC_type(isdw:iedw,jsdw-1:jedw), source=0)
4250
4251 ! Determine the maximum and minimum index range for various directions of OBC points on this PE
4252 ! by first setting these one point outside of the wrong side of the domain.
4253 bt_obc%Is_u_W_obc = iedw + 1 ; bt_obc%Ie_u_W_obc = isdw - 2
4254 bt_obc%js_u_W_obc = jedw + 1 ; bt_obc%je_u_W_obc = jsdw - 1
4255 bt_obc%Is_u_E_obc = iedw + 1 ; bt_obc%Ie_u_E_obc = isdw - 2
4256 bt_obc%js_u_E_obc = jedw + 1 ; bt_obc%je_u_E_obc = jsdw - 1
4257 bt_obc%is_v_S_obc = iedw + 1 ; bt_obc%ie_v_S_obc = isdw - 1
4258 bt_obc%Js_v_S_obc = jedw + 1 ; bt_obc%Je_v_S_obc = jsdw - 2
4259 bt_obc%is_v_N_obc = iedw + 1 ; bt_obc%ie_v_N_obc = isdw - 1
4260 bt_obc%Js_v_N_obc = jedw + 1 ; bt_obc%Je_v_N_obc = jsdw - 2
4261
4262 flather_obc_in_halo = 0
4263 do j=jsdw,jedw ; do i=isdw-1,iedw
4264 bt_obc%u_OBC_type(i,j) = nint(u_obc(i,j))
4265 if (bt_obc%u_OBC_type(i,j) < 0) then ! This point has OBC_DIRECTION_W.
4266 if ((bt_obc%u_OBC_type(i,j) == -flather_obc) .and. (i >= iedw-1)) then
4267 ! There is no need to specify the OBC at this point, but the stencil might need to be increased.
4268 flather_obc_in_halo = 1
4269 else
4270 bt_obc%Is_u_W_obc = min(i, bt_obc%Is_u_W_obc) ; bt_obc%Ie_u_W_obc = max(i, bt_obc%Ie_u_W_obc)
4271 bt_obc%js_u_W_obc = min(j, bt_obc%js_u_W_obc) ; bt_obc%je_u_W_obc = max(j, bt_obc%je_u_W_obc)
4272 endif
4273 endif
4274 if (bt_obc%u_OBC_type(i,j) > 0) then ! This point has OBC_DIRECTION_E.
4275 if ((bt_obc%u_OBC_type(i,j) == flather_obc) .and. (i <= isdw)) then
4276 ! There is no need to specify the OBC at this point, but the stencil might need to be increased.
4277 flather_obc_in_halo = 1
4278 else
4279 bt_obc%Is_u_E_obc = min(i, bt_obc%Is_u_E_obc) ; bt_obc%Ie_u_E_obc = max(i, bt_obc%Ie_u_E_obc)
4280 bt_obc%js_u_E_obc = min(j, bt_obc%js_u_E_obc) ; bt_obc%je_u_E_obc = max(j, bt_obc%je_u_E_obc)
4281 endif
4282 endif
4283 enddo ; enddo
4284
4285 do j=jsdw-1,jedw ; do i=isdw,iedw
4286 bt_obc%v_OBC_type(i,j) = nint(v_obc(i,j))
4287 if (bt_obc%v_OBC_type(i,j) < 0) then ! This point has OBC_DIRECTION_S.
4288 if ((bt_obc%v_OBC_type(i,j) == -flather_obc) .and. (j >= jedw-1)) then
4289 ! There is no need to specify the OBC at this point, but the stencil might need to be increased.
4290 flather_obc_in_halo = 1
4291 else
4292 bt_obc%is_v_S_obc = min(i, bt_obc%is_v_S_obc) ; bt_obc%ie_v_S_obc = max(i, bt_obc%ie_v_S_obc)
4293 bt_obc%Js_v_S_obc = min(j, bt_obc%Js_v_S_obc) ; bt_obc%Je_v_S_obc = max(j, bt_obc%Je_v_S_obc)
4294 endif
4295 endif
4296 if (bt_obc%v_OBC_type(i,j) > 0) then ! This point has OBC_DIRECTION_N.
4297 if ((bt_obc%v_OBC_type(i,j) == flather_obc) .and. (j <= jsdw)) then
4298 ! There is no need to specify the OBC at this point, but the stencil might need to be increased.
4299 flather_obc_in_halo = 1
4300 else
4301 bt_obc%is_v_N_obc = min(i, bt_obc%is_v_N_obc) ; bt_obc%ie_v_N_obc = max(i, bt_obc%ie_v_N_obc)
4302 bt_obc%Js_v_N_obc = min(j, bt_obc%Js_v_N_obc) ; bt_obc%Je_v_N_obc = max(j, bt_obc%Je_v_N_obc)
4303 endif
4304 endif
4305 enddo ; enddo
4306
4307 bt_obc%u_OBCs_on_PE = ((bt_obc%Is_u_E_obc <= iedw) .or. (bt_obc%Is_u_W_obc <= iedw))
4308 bt_obc%v_OBCs_on_PE = ((bt_obc%is_v_N_obc <= iedw) .or. (bt_obc%is_v_S_obc <= iedw))
4309
4310
4311 ! Determine whether there are any OBCs in the opposite halo on any processors in the domain, e.g.,
4312 ! points with OBC_DIRECTION_S in a northern halo.
4313 reversed_obcs = (bt_obc%u_OBCs_on_PE .and. ((bt_obc%Is_u_E_obc <= g%isc-1) .or. (bt_obc%Ie_u_W_obc >= g%iec))) .or. &
4314 (bt_obc%v_OBCs_on_PE .and. ((bt_obc%Js_v_N_obc <= g%jsc-1) .or. (bt_obc%Je_v_S_obc >= g%jec)))
4315 any_reversed_obcs = any_across_pes(reversed_obcs)
4316 if (any_reversed_obcs) call mom_mesg("OBCs in an opposite halo require the use of a wider stencil.", 5)
4317 if (any_reversed_obcs) cs%min_stencil = max(cs%min_stencil, 2)
4318
4319 ! Allocate time-varying arrays that will be used for open boundary conditions.
4320
4321 ! This pair is used with either Flather or specified OBCs.
4322 allocate(bt_obc%ubt_outer(isdw-1:iedw,jsdw:jedw), source=0.0)
4323 allocate(bt_obc%vbt_outer(isdw:iedw,jsdw-1:jedw), source=0.0)
4324 call create_group_pass(bt_obc%pass_uv, bt_obc%ubt_outer, bt_obc%vbt_outer, cs%BT_Domain)
4325
4326 ! This pair is only used with specified OBCs.
4327 allocate(bt_obc%uhbt(isdw-1:iedw,jsdw:jedw), source=0.0)
4328 allocate(bt_obc%vhbt(isdw:iedw,jsdw-1:jedw), source=0.0)
4329 call create_group_pass(bt_obc%pass_uv, bt_obc%uhbt, bt_obc%vhbt, cs%BT_Domain)
4330
4331 if (obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally) then
4332 ! These 3 pairs are only used with Flather OBCs.
4333 allocate(bt_obc%Cg_u(isdw-1:iedw,jsdw:jedw), source=0.0)
4334 allocate(bt_obc%dZ_u(isdw-1:iedw,jsdw:jedw), source=0.0)
4335 allocate(bt_obc%SSH_outer_u(isdw-1:iedw,jsdw:jedw), source=0.0)
4336
4337 allocate(bt_obc%Cg_v(isdw:iedw,jsdw-1:jedw), source=0.0)
4338 allocate(bt_obc%dZ_v(isdw:iedw,jsdw-1:jedw), source=0.0)
4339 allocate(bt_obc%SSH_outer_v(isdw:iedw,jsdw-1:jedw), source=0.0)
4340
4341 call create_group_pass(bt_obc%scalar_pass, bt_obc%SSH_outer_u, bt_obc%SSH_outer_v, cs%BT_Domain, to_all+scalar_pair)
4342 call create_group_pass(bt_obc%scalar_pass, bt_obc%dZ_u, bt_obc%dZ_v, cs%BT_Domain, to_all+scalar_pair)
4343 call create_group_pass(bt_obc%scalar_pass, bt_obc%Cg_u, bt_obc%Cg_v, cs%BT_Domain, to_all+scalar_pair)
4344 endif
4345
4346end subroutine initialize_bt_obc
4347
4348!> This subroutine sets up the time-varying fields in the private structure used to apply the open
4349!! boundary conditions, as developed by Mehmet Ilicak.
4350subroutine set_up_bt_obc(OBC, eta, SpV_avg, BT_OBC, BT_Domain, G, GV, US, CS, MS, halo, use_BT_cont, &
4351 integral_BT_cont, dt_baroclinic, Datu, Datv, BTCL_u, BTCL_v, dgeo_de)
4352 type(ocean_obc_type), target, intent(inout) :: OBC !< An associated pointer to an OBC type.
4353 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of the
4354 !! argument arrays.
4355 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: eta !< The barotropic free surface height anomaly or
4356 !! column mass anomaly [H ~> m or kg m-2].
4357 real, dimension(SZIW_(MS),SZJW_(MS)), intent(in) :: SpV_avg !< The column average specific volume [R-1 ~> m3 kg-1]
4358 type(bt_obc_type), intent(inout) :: BT_OBC !< A structure with the private barotropic arrays
4359 !! related to the open boundary conditions,
4360 !! set by set_up_BT_OBC.
4361 type(mom_domain_type), intent(inout) :: BT_Domain !< MOM_domain_type associated with wide arrays
4362 type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
4363 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
4364 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
4365 type(barotropic_cs), intent(inout) :: CS !< Barotropic control structure
4366 integer, intent(in) :: halo !< The extra halo size to use here.
4367 logical, intent(in) :: use_BT_cont !< If true, use the BT_cont_types to calculate
4368 !! transports.
4369 logical, intent(in) :: integral_BT_cont !< If true, update the barotropic continuity
4370 !! equation directly from the initial condition
4371 !! using the time-integrated barotropic velocity.
4372 real, intent(in) :: dt_baroclinic !< The baroclinic timestep for this cycle of
4373 !! updates to the barotropic solver [T ~> s]
4374 real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: Datu !< A fixed estimate of the face areas at u points
4375 !! [H L ~> m2 or kg m-1].
4376 real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: Datv !< A fixed estimate of the face areas at v points
4377 !! [H L ~> m2 or kg m-1].
4378 type(local_bt_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: BTCL_u !< Structure of information used
4379 !! for a dynamic estimate of the face areas at
4380 !! u-points.
4381 type(local_bt_cont_v_type), dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: BTCL_v !< Structure of information used
4382 !! for a dynamic estimate of the face areas at
4383 !! v-points.
4384 real, intent(in) :: dgeo_de !< The constant of proportionality between
4385 !! geopotential and sea surface height [nondim].
4386 ! Local variables
4387 real :: I_dt ! The inverse of the time interval of this call [T-1 ~> s-1].
4388 integer :: i, j, k, is, ie, js, je, n, nz
4389 integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
4390 integer :: isdw, iedw, jsdw, jedw
4391 type(obc_segment_type), pointer :: segment !< Open boundary segment
4392
4393 is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
4394 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = gv%ke
4395 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
4396 isdw = ms%isdw ; iedw = ms%iedw ; jsdw = ms%jsdw ; jedw = ms%jedw
4397
4398 i_dt = 1.0 / dt_baroclinic
4399
4400 if (bt_obc%u_OBCs_on_PE) then
4401 if (obc%specified_u_BCs_exist_globally) then
4402 do n = 1, obc%number_of_segments
4403 segment => obc%segment(n)
4404 if (segment%is_E_or_W .and. segment%specified) then
4405 do j=segment%HI%jsd,segment%HI%jed ; do i=segment%HI%IsdB,segment%HI%IedB
4406 bt_obc%uhbt(i,j) = 0.
4407 enddo ; enddo
4408 do k=1,nz ; do j=segment%HI%jsd,segment%HI%jed ; do i=segment%HI%IsdB,segment%HI%IedB
4409 bt_obc%uhbt(i,j) = bt_obc%uhbt(i,j) + segment%normal_trans(i,j,k)
4410 enddo ; enddo ; enddo
4411 endif
4412 enddo
4413 endif
4414 do j=js,je ; do i=is-1,ie ; if (bt_obc%u_OBC_type(i,j) /= 0) then
4415 if (abs(bt_obc%u_OBC_type(i,j)) == specified_obc) then ! Eastern or western specified OBC
4416 if (integral_bt_cont) then
4417 bt_obc%ubt_outer(i,j) = uhbt_to_ubt(bt_obc%uhbt(i,j)*dt_baroclinic, btcl_u(i,j)) * i_dt
4418 elseif (use_bt_cont) then
4419 bt_obc%ubt_outer(i,j) = uhbt_to_ubt(bt_obc%uhbt(i,j), btcl_u(i,j))
4420 else
4421 if (datu(i,j) > 0.0) bt_obc%ubt_outer(i,j) = bt_obc%uhbt(i,j) / datu(i,j)
4422 endif
4423 elseif (bt_obc%u_OBC_type(i,j) == flather_obc) then ! Eastern Flather OBC
4424 if (gv%Boussinesq) then
4425 bt_obc%dZ_u(i,j) = cs%bathyT(i,j) + gv%H_to_Z*eta(i,j)
4426 else
4427 bt_obc%dZ_u(i,j) = gv%H_to_RZ * eta(i,j) * spv_avg(i,j)
4428 endif
4429 bt_obc%Cg_u(i,j) = sqrt(dgeo_de * gv%g_prime(1) * bt_obc%dZ_u(i,j))
4430 elseif (bt_obc%u_OBC_type(i,j) == -flather_obc) then ! Western Flather OBC
4431 if (gv%Boussinesq) then
4432 bt_obc%dZ_u(i,j) = cs%bathyT(i+1,j) + gv%H_to_Z*eta(i+1,j)
4433 else
4434 bt_obc%dZ_u(i,j) = gv%H_to_RZ * eta(i+1,j) * spv_avg(i+1,j)
4435 endif
4436 bt_obc%Cg_u(i,j) = sqrt(dgeo_de * gv%g_prime(1) * bt_obc%dZ_u(i,j))
4437 endif
4438 endif ; enddo ; enddo
4439
4440 if (obc%Flather_u_BCs_exist_globally) then
4441 do n = 1, obc%number_of_segments
4442 segment => obc%segment(n)
4443 if (segment%is_E_or_W .and. segment%Flather) then
4444 do j=segment%HI%jsd,segment%HI%jed ; do i=segment%HI%IsdB,segment%HI%IedB
4445 bt_obc%ubt_outer(i,j) = segment%normal_vel_bt(i,j)
4446 bt_obc%SSH_outer_u(i,j) = segment%SSH(i,j) + g%Z_ref
4447 enddo ; enddo
4448 endif
4449 enddo
4450 endif
4451 endif
4452
4453 if (bt_obc%v_OBCs_on_PE) then
4454 if (obc%specified_v_BCs_exist_globally) then
4455 do n = 1, obc%number_of_segments
4456 segment => obc%segment(n)
4457 if (segment%is_N_or_S .and. segment%specified) then
4458 do j=segment%HI%JsdB,segment%HI%JedB ; do i=segment%HI%isd,segment%HI%ied
4459 bt_obc%vhbt(i,j) = 0.
4460 enddo ; enddo
4461 do k=1,nz ; do j=segment%HI%JsdB,segment%HI%JedB ; do i=segment%HI%isd,segment%HI%ied
4462 bt_obc%vhbt(i,j) = bt_obc%vhbt(i,j) + segment%normal_trans(i,j,k)
4463 enddo ; enddo ; enddo
4464 endif
4465 enddo
4466 endif
4467 do j=js-1,je ; do i=is,ie ; if (bt_obc%v_OBC_type(i,j) /= 0) then
4468 if (abs(bt_obc%v_OBC_type(i,j)) == specified_obc) then ! Northern or southern specified OBC
4469 if (integral_bt_cont) then
4470 bt_obc%vbt_outer(i,j) = vhbt_to_vbt(bt_obc%vhbt(i,j)*dt_baroclinic, btcl_v(i,j)) * i_dt
4471 elseif (use_bt_cont) then
4472 bt_obc%vbt_outer(i,j) = vhbt_to_vbt(bt_obc%vhbt(i,j), btcl_v(i,j))
4473 else
4474 if (datv(i,j) > 0.0) bt_obc%vbt_outer(i,j) = bt_obc%vhbt(i,j) / datv(i,j)
4475 endif
4476 elseif (bt_obc%v_OBC_type(i,j) == flather_obc) then ! Northern Flather OBC
4477 if (gv%Boussinesq) then
4478 bt_obc%dZ_v(i,j) = cs%bathyT(i,j) + gv%H_to_Z*eta(i,j)
4479 else
4480 bt_obc%dZ_v(i,j) = gv%H_to_RZ * eta(i,j) * spv_avg(i,j)
4481 endif
4482 bt_obc%Cg_v(i,j) = sqrt(dgeo_de * gv%g_prime(1) * bt_obc%dZ_v(i,j))
4483 elseif (bt_obc%v_OBC_type(i,j) == -flather_obc) then ! Southern Flather OBC
4484 if (gv%Boussinesq) then
4485 bt_obc%dZ_v(i,j) = cs%bathyT(i,j+1) + gv%H_to_Z*eta(i,j+1)
4486 else
4487 bt_obc%dZ_v(i,j) = gv%H_to_RZ * eta(i,j+1) * spv_avg(i,j+1)
4488 endif
4489 bt_obc%Cg_v(i,j) = sqrt(dgeo_de * gv%g_prime(1) * bt_obc%dZ_v(i,j))
4490 endif
4491 endif ; enddo ; enddo
4492 if (obc%Flather_v_BCs_exist_globally) then
4493 do n = 1, obc%number_of_segments
4494 segment => obc%segment(n)
4495 if (segment%is_N_or_S .and. segment%Flather) then
4496 do j=segment%HI%JsdB,segment%HI%JedB ; do i=segment%HI%isd,segment%HI%ied
4497 bt_obc%vbt_outer(i,j) = segment%normal_vel_bt(i,j)
4498 bt_obc%SSH_outer_v(i,j) = segment%SSH(i,j) + g%Z_ref
4499 enddo ; enddo
4500 endif
4501 enddo
4502 endif
4503 endif
4504
4505 call do_group_pass(bt_obc%pass_uv, bt_domain)
4506 if (obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally) &
4507 call do_group_pass(bt_obc%scalar_pass, bt_domain)
4508
4509end subroutine set_up_bt_obc
4510
4511!> Clean up the BT_OBC memory.
4512subroutine destroy_bt_obc(BT_OBC)
4513 type(bt_obc_type), intent(inout) :: BT_OBC !< A structure with the private barotropic arrays
4514 !! related to the open boundary conditions,
4515 !! set by set_up_BT_OBC.
4516
4517 if (allocated(bt_obc%u_OBC_type)) deallocate(bt_obc%u_OBC_type)
4518 if (allocated(bt_obc%v_OBC_type)) deallocate(bt_obc%v_OBC_type)
4519
4520 if (allocated(bt_obc%Cg_u)) deallocate(bt_obc%Cg_u)
4521 if (allocated(bt_obc%dZ_u)) deallocate(bt_obc%dZ_u)
4522 if (allocated(bt_obc%uhbt)) deallocate(bt_obc%uhbt)
4523 if (allocated(bt_obc%ubt_outer)) deallocate(bt_obc%ubt_outer)
4524 if (allocated(bt_obc%SSH_outer_u)) deallocate(bt_obc%SSH_outer_u)
4525
4526 if (allocated(bt_obc%Cg_v)) deallocate(bt_obc%Cg_v)
4527 if (allocated(bt_obc%dZ_v)) deallocate(bt_obc%dZ_v)
4528 if (allocated(bt_obc%vhbt)) deallocate(bt_obc%vhbt)
4529 if (allocated(bt_obc%vbt_outer)) deallocate(bt_obc%vbt_outer)
4530 if (allocated(bt_obc%SSH_outer_v)) deallocate(bt_obc%SSH_outer_v)
4531
4532end subroutine destroy_bt_obc
4533
4534!> btcalc determines the fraction of the total water column in each
4535!! layer at velocity points.
4536subroutine btcalc(h, G, GV, CS, h_u, h_v, may_use_default, OBC)
4537 type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
4538 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
4539 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
4540 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
4541 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
4542 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
4543 optional, intent(in) :: h_u !< The specified effective thicknesses at u-points,
4544 !! perhaps scaled down to account for viscosity and
4545 !! fractional open areas [H ~> m or kg m-2]. These
4546 !! are used here as non-normalized weights for each
4547 !! layer that are converted the normalized weights
4548 !! for determining the barotropic accelerations.
4549 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
4550 optional, intent(in) :: h_v !< The specified effective thicknesses at v-points,
4551 !! perhaps scaled down to account for viscosity and
4552 !! fractional open areas [H ~> m or kg m-2]. These
4553 !! are used here as non-normalized weights for each
4554 !! layer that are converted the normalized weights
4555 !! for determining the barotropic accelerations.
4556 logical, optional, intent(in) :: may_use_default !< An optional logical argument
4557 !! to indicate that the default velocity point
4558 !! thicknesses may be used for this particular
4559 !! calculation, even though the setting of
4560 !! CS%hvel_scheme would usually require that h_u
4561 !! and h_v be passed in.
4562 type(ocean_obc_type), optional, pointer :: obc !< Open boundary control structure.
4563
4564 ! Local variables
4565 real :: hatu(szib_(g),szk_(gv)) ! The layer thicknesses interpolated to u points [H ~> m or kg m-2]
4566 real :: hatv(szi_(g),szk_(gv)) ! The layer thicknesses interpolated to v points [H ~> m or kg m-2]
4567 real :: hatutot(szib_(g)) ! The sum of the layer thicknesses interpolated to u points [H ~> m or kg m-2].
4568 real :: hatvtot(szi_(g)) ! The sum of the layer thicknesses interpolated to v points [H ~> m or kg m-2].
4569 real :: ihatutot(szib_(g)) ! Ihatutot is the inverse of hatutot [H-1 ~> m-1 or m2 kg-1].
4570 real :: ihatvtot(szi_(g)) ! Ihatvtot is the inverse of hatvtot [H-1 ~> m-1 or m2 kg-1].
4571 real :: h_arith ! The arithmetic mean thickness [H ~> m or kg m-2].
4572 real :: h_harm ! The harmonic mean thicknesses [H ~> m or kg m-2].
4573 real :: h_neglect ! A thickness that is so small it is usually lost
4574 ! in roundoff and can be neglected [H ~> m or kg m-2].
4575 real :: wt_arith ! The weight for the arithmetic mean thickness [nondim].
4576 ! The harmonic mean uses a weight of (1 - wt_arith).
4577 real :: e_u(szib_(g),szk_(gv)+1) ! The interface heights at u-velocity points [H ~> m or kg m-2]
4578 real :: e_v(szi_(g),szk_(gv)+1) ! The interface heights at v-velocity points [H ~> m or kg m-2]
4579 real :: d_shallow_u(szi_(g)) ! The height of the shallower of the adjacent bathymetric depths
4580 ! around a u-point (positive upward) [H ~> m or kg m-2]
4581 real :: d_shallow_v(szib_(g))! The height of the shallower of the adjacent bathymetric depths
4582 ! around a v-point (positive upward) [H ~> m or kg m-2]
4583 real :: z_to_h ! A local conversion factor [H Z-1 ~> nondim or kg m-3]
4584
4585 logical :: use_default, test_dflt
4586 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz, i, j, k
4587
4588 if (.not.cs%module_is_initialized) call mom_error(fatal, &
4589 "btcalc: Module MOM_barotropic must be initialized before it is used.")
4590
4591 if (.not.cs%split) return
4592
4593 use_default = .false.
4594 test_dflt = .false. ; if (present(may_use_default)) test_dflt = may_use_default
4595
4596 if (test_dflt) then
4597 if (.not.((present(h_u) .and. present(h_v)) .or. &
4598 (cs%hvel_scheme == harmonic) .or. (cs%hvel_scheme == hybrid) .or.&
4599 (cs%hvel_scheme == arithmetic))) use_default = .true.
4600 else
4601 if (.not.((present(h_u) .and. present(h_v)) .or. &
4602 (cs%hvel_scheme == harmonic) .or. (cs%hvel_scheme == hybrid) .or.&
4603 (cs%hvel_scheme == arithmetic))) call mom_error(fatal, &
4604 "btcalc: Inconsistent settings of optional arguments and hvel_scheme.")
4605 endif
4606
4607 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
4608 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
4609 h_neglect = gv%H_subroundoff
4610
4611
4612 !$OMP parallel do default(none) shared(is,ie,js,je,nz,h_u,CS,h_neglect,h,use_default,G,GV) &
4613 !$OMP private(hatu,hatutot,Ihatutot,e_u,D_shallow_u,h_arith,h_harm,wt_arith,Z_to_H)
4614 do j=js,je
4615 do i=is-1,ie ; hatutot(i) = 0.0 ; enddo
4616
4617 if (present(h_u)) then
4618 do k=1,nz ; do i=is-1,ie
4619 hatu(i,k) = h_u(i,j,k)
4620 hatutot(i) = hatutot(i) + hatu(i,k)
4621 enddo ; enddo
4622 elseif (cs%hvel_scheme == arithmetic) then
4623 do k=1,nz ; do i=is-1,ie
4624 hatu(i,k) = 0.5 * (h(i+1,j,k) + h(i,j,k))
4625 hatutot(i) = hatutot(i) + hatu(i,k)
4626 enddo ; enddo
4627 elseif (cs%hvel_scheme == hybrid .or. use_default) then
4628 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
4629 do i=is-1,ie
4630 e_u(i,nz+1) = -0.5 * z_to_h * (g%bathyT(i+1,j) + g%bathyT(i,j))
4631 d_shallow_u(i) = -z_to_h * min(g%bathyT(i+1,j), g%bathyT(i,j))
4632 enddo
4633 do k=nz,1,-1 ; do i=is-1,ie
4634 e_u(i,k) = e_u(i,k+1) + 0.5 * (h(i+1,j,k) + h(i,j,k))
4635 h_arith = 0.5 * (h(i+1,j,k) + h(i,j,k))
4636 if (e_u(i,k+1) >= d_shallow_u(i)) then
4637 hatu(i,k) = h_arith
4638 else
4639 h_harm = (h(i+1,j,k) * h(i,j,k)) / (h_arith + h_neglect)
4640 if (e_u(i,k) <= d_shallow_u(i)) then
4641 hatu(i,k) = h_harm
4642 else
4643 wt_arith = (e_u(i,k) - d_shallow_u(i)) / (h_arith + h_neglect)
4644 hatu(i,k) = wt_arith*h_arith + (1.0-wt_arith)*h_harm
4645 endif
4646 endif
4647 hatutot(i) = hatutot(i) + hatu(i,k)
4648 enddo ; enddo
4649 elseif (cs%hvel_scheme == harmonic) then
4650 ! Interpolates thicknesses onto u grid points with the
4651 ! second order accurate estimate h = 2*(h+ * h-)/(h+ + h-).
4652 do k=1,nz ; do i=is-1,ie
4653 hatu(i,k) = 2.0*(h(i+1,j,k) * h(i,j,k)) / &
4654 ((h(i+1,j,k) + h(i,j,k)) + h_neglect)
4655 hatutot(i) = hatutot(i) + hatu(i,k)
4656 enddo ; enddo
4657 endif
4658
4659 if (cs%BT_OBC%u_OBCs_on_PE) then
4660 ! Reset velocity point thicknesses and their sums at OBC points
4661 if ((j >= cs%BT_OBC%js_u_E_obc) .and. (j <= cs%BT_OBC%je_u_E_obc)) then
4662 do i = max(is-1,cs%BT_OBC%Is_u_E_obc), min(ie,cs%BT_OBC%Ie_u_E_obc)
4663 if (cs%BT_OBC%u_OBC_type(i,j) > 0) then ! Eastern boundary condition
4664 hatutot(i) = 0.0
4665 do k=1,nz
4666 hatu(i,k) = h(i,j,k)
4667 hatutot(i) = hatutot(i) + hatu(i,k)
4668 enddo
4669 endif
4670 enddo
4671 endif
4672 if ((j >= cs%BT_OBC%js_u_W_obc) .and. (j <= cs%BT_OBC%je_u_W_obc)) then
4673 do i = max(is-1,cs%BT_OBC%Is_u_W_obc), min(ie,cs%BT_OBC%Ie_u_W_obc)
4674 if (cs%BT_OBC%u_OBC_type(i,j) < 0) then ! Western boundary condition
4675 hatutot(i) = 0.0
4676 do k=1,nz
4677 hatu(i,k) = h(i+1,j,k)
4678 hatutot(i) = hatutot(i) + hatu(i,k)
4679 enddo
4680 endif
4681 enddo
4682 endif
4683 endif
4684
4685 ! Determine the fractional thickness of each layer at the velocity points.
4686 do i=is-1,ie ; ihatutot(i) = g%mask2dCu(i,j) / (hatutot(i) + h_neglect) ; enddo
4687 do k=1,nz ; do i=is-1,ie
4688 cs%frhatu(i,j,k) = hatu(i,k) * ihatutot(i)
4689 enddo ; enddo
4690 enddo
4691
4692 !$OMP parallel do default(none) shared(is,ie,js,je,nz,CS,G,GV,h_v,h_neglect,h,use_default) &
4693 !$OMP private(hatv,hatvtot,Ihatvtot,e_v,D_shallow_v,h_arith,h_harm,wt_arith,Z_to_H)
4694 do j=js-1,je
4695 do i=is,ie ; hatvtot(i) = 0.0 ; enddo
4696 if (present(h_v)) then
4697 do k=1,nz ; do i=is,ie
4698 hatv(i,k) = h_v(i,j,k)
4699 hatvtot(i) = hatvtot(i) + hatv(i,k)
4700 enddo ; enddo
4701 elseif (cs%hvel_scheme == arithmetic) then
4702 do k=1,nz ; do i=is,ie
4703 hatv(i,k) = 0.5 * (h(i,j+1,k) + h(i,j,k))
4704 hatvtot(i) = hatvtot(i) + hatv(i,k)
4705 enddo ; enddo
4706 elseif (cs%hvel_scheme == hybrid .or. use_default) then
4707 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
4708 do i=is,ie
4709 e_v(i,nz+1) = -0.5 * z_to_h * (g%bathyT(i,j+1) + g%bathyT(i,j))
4710 d_shallow_v(i) = -z_to_h * min(g%bathyT(i,j+1), g%bathyT(i,j))
4711 enddo
4712 do k=nz,1,-1 ; do i=is,ie
4713 e_v(i,k) = e_v(i,k+1) + 0.5 * (h(i,j+1,k) + h(i,j,k))
4714 h_arith = 0.5 * (h(i,j+1,k) + h(i,j,k))
4715 if (e_v(i,k+1) >= d_shallow_v(i)) then
4716 hatv(i,k) = h_arith
4717 else
4718 h_harm = (h(i,j+1,k) * h(i,j,k)) / (h_arith + h_neglect)
4719 if (e_v(i,k) <= d_shallow_v(i)) then
4720 hatv(i,k) = h_harm
4721 else
4722 wt_arith = (e_v(i,k) - d_shallow_v(i)) / (h_arith + h_neglect)
4723 hatv(i,k) = wt_arith*h_arith + (1.0-wt_arith)*h_harm
4724 endif
4725 endif
4726 hatvtot(i) = hatvtot(i) + hatv(i,k)
4727 enddo ; enddo
4728 elseif (cs%hvel_scheme == harmonic) then
4729 do k=1,nz ; do i=is,ie
4730 hatv(i,k) = 2.0*(h(i,j+1,k) * h(i,j,k)) / &
4731 ((h(i,j+1,k) + h(i,j,k)) + h_neglect)
4732 hatvtot(i) = hatvtot(i) + hatv(i,k)
4733 enddo ; enddo
4734 endif
4735
4736 if (cs%BT_OBC%v_OBCs_on_PE) then
4737 ! Reset v-velocity point thicknesses and their sums at OBC points
4738 if ((j >= cs%BT_OBC%Js_v_N_obc) .and. (j <= cs%BT_OBC%Je_v_N_obc)) then
4739 do i = max(is,cs%BT_OBC%is_v_N_obc), min(ie,cs%BT_OBC%ie_v_N_obc)
4740 if (cs%BT_OBC%v_OBC_type(i,j) > 0) then ! Northern boundary condition
4741 hatvtot(i) = 0.0
4742 do k=1,nz
4743 hatv(i,k) = h(i,j,k)
4744 hatvtot(i) = hatvtot(i) + hatv(i,k)
4745 enddo
4746 endif
4747 enddo
4748 endif
4749 if ((j >= cs%BT_OBC%Js_v_S_obc) .and. (j <= cs%BT_OBC%Je_v_S_obc)) then
4750 do i = max(is,cs%BT_OBC%is_v_S_obc), min(ie,cs%BT_OBC%ie_v_S_obc)
4751 if (cs%BT_OBC%v_OBC_type(i,j) < 0) then ! Southern boundary condition
4752 hatvtot(i) = 0.0
4753 do k=1,nz
4754 hatv(i,k) = h(i,j+1,k)
4755 hatvtot(i) = hatvtot(i) + hatv(i,k)
4756 enddo
4757 endif
4758 enddo
4759 endif
4760 endif
4761
4762 ! Determine the fractional thickness of each layer at the velocity points.
4763 do i=is,ie ; ihatvtot(i) = g%mask2dCv(i,j) / (hatvtot(i) + h_neglect) ; enddo
4764 do k=1,nz ; do i=is,ie
4765 cs%frhatv(i,j,k) = hatv(i,k) * ihatvtot(i)
4766 enddo ; enddo
4767 enddo
4768
4769 if (cs%debug) then
4770 call uvchksum("btcalc frhat[uv]", cs%frhatu, cs%frhatv, g%HI, &
4771 haloshift=0, symmetric=.true., omit_corners=.true., &
4772 scalar_pair=.true.)
4773 if (present(h_u) .and. present(h_v)) &
4774 call uvchksum("btcalc h_[uv]", h_u, h_v, g%HI, haloshift=0, &
4775 symmetric=.true., omit_corners=.true., unscale=gv%H_to_MKS, &
4776 scalar_pair=.true.)
4777 call hchksum(h, "btcalc h", g%HI, haloshift=1, unscale=gv%H_to_MKS)
4778 endif
4779
4780end subroutine btcalc
4781
4782!> The function find_uhbt determines the zonal transport for a given velocity, or with
4783!! INTEGRAL_BT_CONT=True it determines the time-integrated zonal transport for a given
4784!! time-integrated velocity.
4785function find_uhbt(u, BTC) result(uhbt)
4786 real, intent(in) :: u !< The local zonal velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]
4787 type(local_bt_cont_u_type), intent(in) :: btc !< A structure containing various fields that
4788 !! allow the barotropic transports to be calculated consistently
4789 !! with the layers' continuity equations. The dimensions of some
4790 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4791
4792 real :: uhbt !< The zonal barotropic transport [L2 H T-1 ~> m3 s-1] or time integrated transport [L2 H ~> m3]
4793
4794 if (u == 0.0) then
4795 uhbt = 0.0
4796 elseif (u < btc%uBT_EE) then
4797 uhbt = (u - btc%uBT_EE) * btc%FA_u_EE + btc%uh_EE
4798 elseif (u < 0.0) then
4799 uhbt = u * (btc%FA_u_E0 + btc%uh_crvE * u**2)
4800 elseif (u <= btc%uBT_WW) then
4801 uhbt = u * (btc%FA_u_W0 + btc%uh_crvW * u**2)
4802 else ! (u > BTC%uBT_WW)
4803 uhbt = (u - btc%uBT_WW) * btc%FA_u_WW + btc%uh_WW
4804 endif
4805
4806end function find_uhbt
4807
4808!> The function find_duhbt_du determines the marginal zonal face area for a given velocity, or
4809!! with INTEGRAL_BT_CONT=True for a given time-integrated velocity.
4810function find_duhbt_du(u, BTC) result(duhbt_du)
4811 real, intent(in) :: u !< The local zonal velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]
4812 type(local_bt_cont_u_type), intent(in) :: btc !< A structure containing various fields that
4813 !! allow the barotropic transports to be calculated consistently
4814 !! with the layers' continuity equations. The dimensions of some
4815 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4816 real :: duhbt_du !< The zonal barotropic face area [L H ~> m2 or kg m-1]
4817
4818 if (u == 0.0) then
4819 duhbt_du = 0.5*(btc%FA_u_E0 + btc%FA_u_W0) ! Note the potential discontinuity here.
4820 elseif (u < btc%uBT_EE) then
4821 duhbt_du = btc%FA_u_EE
4822 elseif (u < 0.0) then
4823 duhbt_du = (btc%FA_u_E0 + 3.0*btc%uh_crvE * u**2)
4824 elseif (u <= btc%uBT_WW) then
4825 duhbt_du = (btc%FA_u_W0 + 3.0*btc%uh_crvW * u**2)
4826 else ! (u > BTC%uBT_WW)
4827 duhbt_du = btc%FA_u_WW
4828 endif
4829
4830end function find_duhbt_du
4831
4832!> This function inverts the transport function to determine the barotopic
4833!! velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True
4834!! this finds the time-integrated velocity that is consistent with a time-integrated transport.
4835function uhbt_to_ubt(uhbt, BTC) result(ubt)
4836 real, intent(in) :: uhbt !< The barotropic zonal transport that should be inverted for,
4837 !! [H L2 T-1 ~> m3 s-1 or kg s-1] or the time-integrated
4838 !! transport [H L2 ~> m3 or kg].
4839 type(local_bt_cont_u_type), intent(in) :: btc !< A structure containing various fields that allow the
4840 !! barotropic transports to be calculated consistently with the
4841 !! layers' continuity equations. The dimensions of some
4842 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4843 real :: ubt !< The result - The velocity that gives uhbt transport [L T-1 ~> m s-1]
4844 !! or the time-integrated velocity [L ~> m].
4845
4846 ! Local variables
4847 real :: ubt_min, ubt_max ! Bounding values of vbt [L T-1 ~> m s-1] or [L ~> m]
4848 real :: uhbt_err ! The transport error [H L2 T-1 ~> m3 s-1 or kg s-1] or [H L2 ~> m3 or kg].
4849 real :: derr_du ! The change in transport error with vbt, i.e. the face area [H L ~> m2 or kg m-1].
4850 real :: uherr_min, uherr_max ! The bounding values of the transport error [H L2 T-1 ~> m3 s-1 or kg s-1]
4851 ! or [H L2 ~> m3 or kg].
4852 real, parameter :: tol = 1.0e-10 ! A fractional match tolerance [nondim]
4853 real, parameter :: vs1 = 1.25 ! Nondimensional parameters used in limiting
4854 real, parameter :: vs2 = 2.0 ! the velocity, starting at vs1, with the
4855 ! maximum increase of vs2, both [nondim].
4856 integer :: itt, max_itt = 20
4857
4858 ! Find the value of ubt that gives uhbt.
4859 if (uhbt == 0.0) then
4860 ubt = 0.0
4861 elseif (uhbt < btc%uh_EE) then
4862 ubt = btc%uBT_EE + (uhbt - btc%uh_EE) / btc%FA_u_EE
4863 elseif (uhbt < 0.0) then
4864 ! Iterate to convergence with Newton's method (when bounded) and the
4865 ! false position method otherwise. ubt will be negative.
4866 ubt_min = btc%uBT_EE ; uherr_min = btc%uh_EE - uhbt
4867 ubt_max = 0.0 ; uherr_max = -uhbt
4868 ! Use a false-position method first guess.
4869 ubt = btc%uBT_EE * (uhbt / btc%uh_EE)
4870 do itt = 1, max_itt
4871 uhbt_err = ubt * (btc%FA_u_E0 + btc%uh_crvE * ubt**2) - uhbt
4872
4873 if (abs(uhbt_err) < tol*abs(uhbt)) exit
4874 if (uhbt_err > 0.0) then ; ubt_max = ubt ; uherr_max = uhbt_err ; endif
4875 if (uhbt_err < 0.0) then ; ubt_min = ubt ; uherr_min = uhbt_err ; endif
4876
4877 derr_du = btc%FA_u_E0 + 3.0 * btc%uh_crvE * ubt**2
4878 if ((uhbt_err >= derr_du*(ubt - ubt_min)) .or. &
4879 (-uhbt_err >= derr_du*(ubt_max - ubt)) .or. (derr_du <= 0.0)) then
4880 ! Use a false-position method guess.
4881 ubt = ubt_max + (ubt_min-ubt_max) * (uherr_max / (uherr_max-uherr_min))
4882 else ! Use Newton's method.
4883 ubt = ubt - uhbt_err / derr_du
4884 if (abs(uhbt_err) < (0.01*tol)*abs(ubt_min*derr_du)) exit
4885 endif
4886 enddo
4887 elseif (uhbt <= btc%uh_WW) then
4888 ! Iterate to convergence with Newton's method. ubt will be positive.
4889 ubt_min = 0.0 ; uherr_min = -uhbt
4890 ubt_max = btc%uBT_WW ; uherr_max = btc%uh_WW - uhbt
4891 ! Use a false-position method first guess.
4892 ubt = btc%uBT_WW * (uhbt / btc%uh_WW)
4893 do itt = 1, max_itt
4894 uhbt_err = ubt * (btc%FA_u_W0 + btc%uh_crvW * ubt**2) - uhbt
4895
4896 if (abs(uhbt_err) < tol*abs(uhbt)) exit
4897 if (uhbt_err > 0.0) then ; ubt_max = ubt ; uherr_max = uhbt_err ; endif
4898 if (uhbt_err < 0.0) then ; ubt_min = ubt ; uherr_min = uhbt_err ; endif
4899
4900 derr_du = btc%FA_u_W0 + 3.0 * btc%uh_crvW * ubt**2
4901 if ((uhbt_err >= derr_du*(ubt - ubt_min)) .or. &
4902 (-uhbt_err >= derr_du*(ubt_max - ubt)) .or. (derr_du <= 0.0)) then
4903 ! Use a false-position method guess.
4904 ubt = ubt_min + (ubt_max-ubt_min) * (-uherr_min / (uherr_max-uherr_min))
4905 else ! Use Newton's method.
4906 ubt = ubt - uhbt_err / derr_du
4907 if (abs(uhbt_err) < (0.01*tol)*(ubt_max*derr_du)) exit
4908 endif
4909 enddo
4910 else ! (uhbt > BTC%uh_WW)
4911 ubt = btc%uBT_WW + (uhbt - btc%uh_WW) / btc%FA_u_WW
4912 endif
4913
4914end function uhbt_to_ubt
4915
4916!> The function find_vhbt determines the meridional transport for a given velocity, or with
4917!! INTEGRAL_BT_CONT=True it determines the time-integrated meridional transport for a given
4918!! time-integrated velocity.
4919function find_vhbt(v, BTC) result(vhbt)
4920 real, intent(in) :: v !< The local meridional velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]
4921 type(local_bt_cont_v_type), intent(in) :: btc !< A structure containing various fields that
4922 !! allow the barotropic transports to be calculated consistently
4923 !! with the layers' continuity equations. The dimensions of some
4924 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4925 real :: vhbt !< The meridional barotropic transport [L2 H T-1 ~> m3 s-1] or time integrated transport [L2 H ~> m3]
4926
4927 if (v == 0.0) then
4928 vhbt = 0.0
4929 elseif (v < btc%vBT_NN) then
4930 vhbt = (v - btc%vBT_NN) * btc%FA_v_NN + btc%vh_NN
4931 elseif (v < 0.0) then
4932 vhbt = v * (btc%FA_v_N0 + btc%vh_crvN * v**2)
4933 elseif (v <= btc%vBT_SS) then
4934 vhbt = v * (btc%FA_v_S0 + btc%vh_crvS * v**2)
4935 else ! (v > BTC%vBT_SS)
4936 vhbt = (v - btc%vBT_SS) * btc%FA_v_SS + btc%vh_SS
4937 endif
4938
4939end function find_vhbt
4940
4941!> The function find_dvhbt_dv determines the marginal meridional face area for a given velocity, or
4942!! with INTEGRAL_BT_CONT=True for a given time-integrated velocity.
4943function find_dvhbt_dv(v, BTC) result(dvhbt_dv)
4944 real, intent(in) :: v !< The local meridional velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]
4945 type(local_bt_cont_v_type), intent(in) :: btc !< A structure containing various fields that
4946 !! allow the barotropic transports to be calculated consistently
4947 !! with the layers' continuity equations. The dimensions of some
4948 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4949 real :: dvhbt_dv !< The meridional barotropic face area [L H ~> m2 or kg m-1]
4950
4951 if (v == 0.0) then
4952 dvhbt_dv = 0.5*(btc%FA_v_N0 + btc%FA_v_S0) ! Note the potential discontinuity here.
4953 elseif (v < btc%vBT_NN) then
4954 dvhbt_dv = btc%FA_v_NN
4955 elseif (v < 0.0) then
4956 dvhbt_dv = btc%FA_v_N0 + 3.0*btc%vh_crvN * v**2
4957 elseif (v <= btc%vBT_SS) then
4958 dvhbt_dv = btc%FA_v_S0 + 3.0*btc%vh_crvS * v**2
4959 else ! (v > BTC%vBT_SS)
4960 dvhbt_dv = btc%FA_v_SS
4961 endif
4962
4963end function find_dvhbt_dv
4964
4965!> This function inverts the transport function to determine the barotopic
4966!! velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True
4967!! this finds the time-integrated velocity that is consistent with a time-integrated transport.
4968function vhbt_to_vbt(vhbt, BTC) result(vbt)
4969 real, intent(in) :: vhbt !< The barotropic meridional transport that should be
4970 !! inverted for [H L2 T-1 ~> m3 s-1 or kg s-1] or the
4971 !! time-integrated transport [H L2 ~> m3 or kg].
4972 type(local_bt_cont_v_type), intent(in) :: btc !< A structure containing various fields that allow the
4973 !! barotropic transports to be calculated consistently
4974 !! with the layers' continuity equations. The dimensions of some
4975 !! of the elements in this type vary depending on INTEGRAL_BT_CONT.
4976 real :: vbt !< The result - The velocity that gives vhbt transport [L T-1 ~> m s-1]
4977 !! or the time-integrated velocity [L ~> m].
4978
4979 ! Local variables
4980 real :: vbt_min, vbt_max ! Bounding values of vbt [L T-1 ~> m s-1] or [L ~> m]
4981 real :: vhbt_err ! The transport error [H L2 T-1 ~> m3 s-1 or kg s-1] or [H L2 ~> m3 or kg].
4982 real :: derr_dv ! The change in transport error with vbt, i.e. the face area [H L ~> m2 or kg m-1].
4983 real :: vherr_min, vherr_max ! The bounding values of the transport error [H L2 T-1 ~> m3 s-1 or kg s-1]
4984 ! or [H L2 ~> m3 or kg].
4985 real, parameter :: tol = 1.0e-10 ! A fractional match tolerance [nondim]
4986 real, parameter :: vs1 = 1.25 ! Nondimensional parameters used in limiting
4987 real, parameter :: vs2 = 2.0 ! the velocity, starting at vs1, with the
4988 ! maximum increase of vs2, both [nondim].
4989 integer :: itt, max_itt = 20
4990
4991 ! Find the value of vbt that gives vhbt.
4992 if (vhbt == 0.0) then
4993 vbt = 0.0
4994 elseif (vhbt < btc%vh_NN) then
4995 vbt = btc%vBT_NN + (vhbt - btc%vh_NN) / btc%FA_v_NN
4996 elseif (vhbt < 0.0) then
4997 ! Iterate to convergence with Newton's method (when bounded) and the
4998 ! false position method otherwise. vbt will be negative.
4999 vbt_min = btc%vBT_NN ; vherr_min = btc%vh_NN - vhbt
5000 vbt_max = 0.0 ; vherr_max = -vhbt
5001 ! Use a false-position method first guess.
5002 vbt = btc%vBT_NN * (vhbt / btc%vh_NN)
5003 do itt = 1, max_itt
5004 vhbt_err = vbt * (btc%FA_v_N0 + btc%vh_crvN * vbt**2) - vhbt
5005
5006 if (abs(vhbt_err) < tol*abs(vhbt)) exit
5007 if (vhbt_err > 0.0) then ; vbt_max = vbt ; vherr_max = vhbt_err ; endif
5008 if (vhbt_err < 0.0) then ; vbt_min = vbt ; vherr_min = vhbt_err ; endif
5009
5010 derr_dv = btc%FA_v_N0 + 3.0 * btc%vh_crvN * vbt**2
5011 if ((vhbt_err >= derr_dv*(vbt - vbt_min)) .or. &
5012 (-vhbt_err >= derr_dv*(vbt_max - vbt)) .or. (derr_dv <= 0.0)) then
5013 ! Use a false-position method guess.
5014 vbt = vbt_max + (vbt_min-vbt_max) * (vherr_max / (vherr_max-vherr_min))
5015 else ! Use Newton's method.
5016 vbt = vbt - vhbt_err / derr_dv
5017 if (abs(vhbt_err) < (0.01*tol)*abs(derr_dv*vbt_min)) exit
5018 endif
5019 enddo
5020 elseif (vhbt <= btc%vh_SS) then
5021 ! Iterate to convergence with Newton's method. vbt will be positive.
5022 vbt_min = 0.0 ; vherr_min = -vhbt
5023 vbt_max = btc%vBT_SS ; vherr_max = btc%vh_SS - vhbt
5024 ! Use a false-position method first guess.
5025 vbt = btc%vBT_SS * (vhbt / btc%vh_SS)
5026 do itt = 1, max_itt
5027 vhbt_err = vbt * (btc%FA_v_S0 + btc%vh_crvS * vbt**2) - vhbt
5028
5029 if (abs(vhbt_err) < tol*abs(vhbt)) exit
5030 if (vhbt_err > 0.0) then ; vbt_max = vbt ; vherr_max = vhbt_err ; endif
5031 if (vhbt_err < 0.0) then ; vbt_min = vbt ; vherr_min = vhbt_err ; endif
5032
5033 derr_dv = btc%FA_v_S0 + 3.0 * btc%vh_crvS * vbt**2
5034 if ((vhbt_err >= derr_dv*(vbt - vbt_min)) .or. &
5035 (-vhbt_err >= derr_dv*(vbt_max - vbt)) .or. (derr_dv <= 0.0)) then
5036 ! Use a false-position method guess.
5037 vbt = vbt_min + (vbt_max-vbt_min) * (-vherr_min / (vherr_max-vherr_min))
5038 else ! Use Newton's method.
5039 vbt = vbt - vhbt_err / derr_dv
5040 if (abs(vhbt_err) < (0.01*tol)*(vbt_max*derr_dv)) exit
5041 endif
5042 enddo
5043 else ! (vhbt > BTC%vh_SS)
5044 vbt = btc%vBT_SS + (vhbt - btc%vh_SS) / btc%FA_v_SS
5045 endif
5046
5047end function vhbt_to_vbt
5048
5049!> This subroutine sets up reordered versions of the BT_cont type in the
5050!! local_BT_cont types, which have wide halos properly filled in.
5051subroutine set_local_bt_cont_types(BT_cont, BTCL_u, BTCL_v, G, US, MS, BT_Domain, halo, dt_baroclinic)
5052 type(bt_cont_type), intent(inout) :: BT_cont !< The BT_cont_type input to the barotropic solver
5053 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of
5054 !! the argument arrays
5055 type(local_bt_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), &
5056 intent(out) :: BTCL_u !< A structure with the u information from BT_cont
5057 type(local_bt_cont_v_type), dimension(SZIW_(MS),SZJBW_(MS)), &
5058 intent(out) :: BTCL_v !< A structure with the v information from BT_cont
5059 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
5060 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
5061 type(mom_domain_type), intent(inout) :: BT_Domain !< The domain to use for updating the halos
5062 !! of wide arrays
5063 integer, intent(in) :: halo !< The extra halo size to use here
5064 real, optional, intent(in) :: dt_baroclinic !< The baroclinic time step [T ~> s], which
5065 !! is provided if INTEGRAL_BT_CONTINUITY is true.
5066
5067 ! Local variables
5068 real, dimension(SZIBW_(MS),SZJW_(MS)) :: &
5069 u_polarity, & ! An array used to test for halo update polarity [nondim]
5070 uBT_EE, uBT_WW, & ! Zonal velocities at which the form of the fit changes [L T-1 ~> m s-1]
5071 FA_u_EE, FA_u_E0, FA_u_W0, FA_u_WW ! Zonal face areas [H L ~> m2 or kg m-1]
5072 real, dimension(SZIW_(MS),SZJBW_(MS)) :: &
5073 v_polarity, & ! An array used to test for halo update polarity [nondim]
5074 vBT_NN, vBT_SS, & ! Meridional velocities at which the form of the fit changes [L T-1 ~> m s-1]
5075 FA_v_NN, FA_v_N0, FA_v_S0, FA_v_SS ! Meridional face areas [H L ~> m2 or kg m-1]
5076 real :: dt ! The baroclinic timestep [T ~> s] or 1.0 [nondim]
5077 real, parameter :: C1_3 = 1.0/3.0 ! [nondim]
5078 integer :: i, j, is, ie, js, je, hs
5079
5080 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
5081 hs = max(halo,0)
5082 dt = 1.0 ; if (present(dt_baroclinic)) dt = dt_baroclinic
5083
5084 ! Copy the BT_cont arrays into symmetric, potentially wide haloed arrays.
5085 !$OMP parallel default(shared)
5086 !$OMP do
5087 do j=js-hs,je+hs ; do i=is-hs-1,ie+hs
5088 u_polarity(i,j) = 1.0
5089 ubt_ee(i,j) = 0.0 ; ubt_ww(i,j) = 0.0
5090 fa_u_ee(i,j) = 0.0 ; fa_u_e0(i,j) = 0.0 ; fa_u_w0(i,j) = 0.0 ; fa_u_ww(i,j) = 0.0
5091 enddo ; enddo
5092 !$OMP do
5093 do j=js-hs-1,je+hs ; do i=is-hs,ie+hs
5094 v_polarity(i,j) = 1.0
5095 vbt_nn(i,j) = 0.0 ; vbt_ss(i,j) = 0.0
5096 fa_v_nn(i,j) = 0.0 ; fa_v_n0(i,j) = 0.0 ; fa_v_s0(i,j) = 0.0 ; fa_v_ss(i,j) = 0.0
5097 enddo ; enddo
5098 !$OMP do
5099 do j=js,je ; do i=is-1,ie
5100 ubt_ee(i,j) = bt_cont%uBT_EE(i,j) ; ubt_ww(i,j) = bt_cont%uBT_WW(i,j)
5101 fa_u_ee(i,j) = bt_cont%FA_u_EE(i,j) ; fa_u_e0(i,j) = bt_cont%FA_u_E0(i,j)
5102 fa_u_w0(i,j) = bt_cont%FA_u_W0(i,j) ; fa_u_ww(i,j) = bt_cont%FA_u_WW(i,j)
5103 enddo ; enddo
5104 !$OMP do
5105 do j=js-1,je ; do i=is,ie
5106 vbt_nn(i,j) = bt_cont%vBT_NN(i,j) ; vbt_ss(i,j) = bt_cont%vBT_SS(i,j)
5107 fa_v_nn(i,j) = bt_cont%FA_v_NN(i,j) ; fa_v_n0(i,j) = bt_cont%FA_v_N0(i,j)
5108 fa_v_s0(i,j) = bt_cont%FA_v_S0(i,j) ; fa_v_ss(i,j) = bt_cont%FA_v_SS(i,j)
5109 enddo ; enddo
5110 !$OMP end parallel
5111
5112 if (id_clock_calc_pre > 0) call cpu_clock_end(id_clock_calc_pre)
5113 if (id_clock_pass_pre > 0) call cpu_clock_begin(id_clock_pass_pre)
5114!--- begin setup for group halo update
5115 call create_group_pass(bt_cont%pass_polarity_BT, u_polarity, v_polarity, bt_domain)
5116 call create_group_pass(bt_cont%pass_polarity_BT, ubt_ee, vbt_nn, bt_domain)
5117 call create_group_pass(bt_cont%pass_polarity_BT, ubt_ww, vbt_ss, bt_domain)
5118
5119 call create_group_pass(bt_cont%pass_FA_uv, fa_u_ee, fa_v_nn, bt_domain, to_all+scalar_pair)
5120 call create_group_pass(bt_cont%pass_FA_uv, fa_u_e0, fa_v_n0, bt_domain, to_all+scalar_pair)
5121 call create_group_pass(bt_cont%pass_FA_uv, fa_u_w0, fa_v_s0, bt_domain, to_all+scalar_pair)
5122 call create_group_pass(bt_cont%pass_FA_uv, fa_u_ww, fa_v_ss, bt_domain, to_all+scalar_pair)
5123!--- end setup for group halo update
5124 ! Do halo updates on BT_cont.
5125 call do_group_pass(bt_cont%pass_polarity_BT, bt_domain)
5126 call do_group_pass(bt_cont%pass_FA_uv, bt_domain)
5127 if (id_clock_pass_pre > 0) call cpu_clock_end(id_clock_pass_pre)
5128 if (id_clock_calc_pre > 0) call cpu_clock_begin(id_clock_calc_pre)
5129
5130 !$OMP parallel default(shared)
5131 !$OMP do
5132 do j=js-hs,je+hs ; do i=is-hs-1,ie+hs
5133 btcl_u(i,j)%FA_u_EE = fa_u_ee(i,j) ; btcl_u(i,j)%FA_u_E0 = fa_u_e0(i,j)
5134 btcl_u(i,j)%FA_u_W0 = fa_u_w0(i,j) ; btcl_u(i,j)%FA_u_WW = fa_u_ww(i,j)
5135 btcl_u(i,j)%uBT_EE = dt*ubt_ee(i,j) ; btcl_u(i,j)%uBT_WW = dt*ubt_ww(i,j)
5136 ! Check for reversed polarity in the tripolar halo regions.
5137 if (u_polarity(i,j) < 0.0) then
5138 call swap(btcl_u(i,j)%FA_u_EE, btcl_u(i,j)%FA_u_WW)
5139 call swap(btcl_u(i,j)%FA_u_E0, btcl_u(i,j)%FA_u_W0)
5140 call swap(btcl_u(i,j)%uBT_EE, btcl_u(i,j)%uBT_WW)
5141 endif
5142
5143 btcl_u(i,j)%uh_EE = btcl_u(i,j)%uBT_EE * &
5144 (c1_3 * (2.0*btcl_u(i,j)%FA_u_E0 + btcl_u(i,j)%FA_u_EE))
5145 btcl_u(i,j)%uh_WW = btcl_u(i,j)%uBT_WW * &
5146 (c1_3 * (2.0*btcl_u(i,j)%FA_u_W0 + btcl_u(i,j)%FA_u_WW))
5147
5148 btcl_u(i,j)%uh_crvE = 0.0 ; btcl_u(i,j)%uh_crvW = 0.0
5149 if (abs(btcl_u(i,j)%uBT_WW) > 0.0) btcl_u(i,j)%uh_crvW = &
5150 (c1_3 * (btcl_u(i,j)%FA_u_WW - btcl_u(i,j)%FA_u_W0)) / btcl_u(i,j)%uBT_WW**2
5151 if (abs(btcl_u(i,j)%uBT_EE) > 0.0) btcl_u(i,j)%uh_crvE = &
5152 (c1_3 * (btcl_u(i,j)%FA_u_EE - btcl_u(i,j)%FA_u_E0)) / btcl_u(i,j)%uBT_EE**2
5153 enddo ; enddo
5154 !$OMP do
5155 do j=js-hs-1,je+hs ; do i=is-hs,ie+hs
5156 btcl_v(i,j)%FA_v_NN = fa_v_nn(i,j) ; btcl_v(i,j)%FA_v_N0 = fa_v_n0(i,j)
5157 btcl_v(i,j)%FA_v_S0 = fa_v_s0(i,j) ; btcl_v(i,j)%FA_v_SS = fa_v_ss(i,j)
5158 btcl_v(i,j)%vBT_NN = dt*vbt_nn(i,j) ; btcl_v(i,j)%vBT_SS = dt*vbt_ss(i,j)
5159 ! Check for reversed polarity in the tripolar halo regions.
5160 if (v_polarity(i,j) < 0.0) then
5161 call swap(btcl_v(i,j)%FA_v_NN, btcl_v(i,j)%FA_v_SS)
5162 call swap(btcl_v(i,j)%FA_v_N0, btcl_v(i,j)%FA_v_S0)
5163 call swap(btcl_v(i,j)%vBT_NN, btcl_v(i,j)%vBT_SS)
5164 endif
5165
5166 btcl_v(i,j)%vh_NN = btcl_v(i,j)%vBT_NN * &
5167 (c1_3 * (2.0*btcl_v(i,j)%FA_v_N0 + btcl_v(i,j)%FA_v_NN))
5168 btcl_v(i,j)%vh_SS = btcl_v(i,j)%vBT_SS * &
5169 (c1_3 * (2.0*btcl_v(i,j)%FA_v_S0 + btcl_v(i,j)%FA_v_SS))
5170
5171 btcl_v(i,j)%vh_crvN = 0.0 ; btcl_v(i,j)%vh_crvS = 0.0
5172 if (abs(btcl_v(i,j)%vBT_SS) > 0.0) btcl_v(i,j)%vh_crvS = &
5173 (c1_3 * (btcl_v(i,j)%FA_v_SS - btcl_v(i,j)%FA_v_S0)) / btcl_v(i,j)%vBT_SS**2
5174 if (abs(btcl_v(i,j)%vBT_NN) > 0.0) btcl_v(i,j)%vh_crvN = &
5175 (c1_3 * (btcl_v(i,j)%FA_v_NN - btcl_v(i,j)%FA_v_N0)) / btcl_v(i,j)%vBT_NN**2
5176 enddo ; enddo
5177 !$OMP end parallel
5178end subroutine set_local_bt_cont_types
5179
5180
5181!> Adjust_local_BT_cont_types expands the range of velocities with a cubic curve
5182!! translating velocities into transports to match the initial values of velocities and
5183!! summed transports when the velocities are larger than the first guesses of the cubic
5184!! transition velocities used to set up the local_BT_cont types.
5185subroutine adjust_local_bt_cont_types(ubt, uhbt, vbt, vhbt, BTCL_u, BTCL_v, &
5186 G, US, MS, halo, dt_baroclinic)
5187 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of the argument arrays.
5188 real, dimension(SZIBW_(MS),SZJW_(MS)), &
5189 intent(in) :: ubt !< The linearization zonal barotropic velocity [L T-1 ~> m s-1].
5190 real, dimension(SZIBW_(MS),SZJW_(MS)), &
5191 intent(in) :: uhbt !< The linearization zonal barotropic transport
5192 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
5193 real, dimension(SZIW_(MS),SZJBW_(MS)), &
5194 intent(in) :: vbt !< The linearization meridional barotropic velocity [L T-1 ~> m s-1].
5195 real, dimension(SZIW_(MS),SZJBW_(MS)), &
5196 intent(in) :: vhbt !< The linearization meridional barotropic transport
5197 !! [H L2 T-1 ~> m3 s-1 or kg s-1].
5198 type(local_bt_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), &
5199 intent(out) :: BTCL_u !< A structure with the u information from BT_cont.
5200 type(local_bt_cont_v_type), dimension(SZIW_(MS),SZJBW_(MS)), &
5201 intent(out) :: BTCL_v !< A structure with the v information from BT_cont.
5202 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
5203 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
5204 integer, intent(in) :: halo !< The extra halo size to use here.
5205 real, optional, intent(in) :: dt_baroclinic !< The baroclinic time step [T ~> s], which is
5206 !! provided if INTEGRAL_BT_CONTINUITY is true.
5207
5208 ! Local variables
5209 real :: dt ! The baroclinic timestep [T ~> s] or 1.0 [nondim]
5210 real, parameter :: C1_3 = 1.0/3.0 ! [nondim]
5211 integer :: i, j, is, ie, js, je, hs
5212
5213 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
5214 hs = max(halo,0)
5215 dt = 1.0 ; if (present(dt_baroclinic)) dt = dt_baroclinic
5216
5217 !$OMP parallel do default(shared)
5218 do j=js-hs,je+hs ; do i=is-hs-1,ie+hs
5219 if ((dt*ubt(i,j) > btcl_u(i,j)%uBT_WW) .and. (dt*uhbt(i,j) > btcl_u(i,j)%uh_WW)) then
5220 ! Expand the cubic fit to use this new point. ubt is negative.
5221 btcl_u(i,j)%ubt_WW = dt * ubt(i,j)
5222 if (3.0*uhbt(i,j) < 2.0*ubt(i,j) * btcl_u(i,j)%FA_u_W0) then
5223 ! No further bounding is needed.
5224 btcl_u(i,j)%uh_crvW = (uhbt(i,j) - ubt(i,j) * btcl_u(i,j)%FA_u_W0) / (dt**2 * ubt(i,j)**3)
5225 else ! This should not happen often!
5226 btcl_u(i,j)%FA_u_W0 = 1.5*uhbt(i,j) / ubt(i,j)
5227 btcl_u(i,j)%uh_crvW = -0.5*uhbt(i,j) / (dt**2 * ubt(i,j)**3)
5228 endif
5229 btcl_u(i,j)%uh_WW = dt * uhbt(i,j)
5230 ! I don't know whether this is helpful.
5231! BTCL_u(I,j)%FA_u_WW = min(BTCL_u(I,j)%FA_u_WW, uhbt(I,j) / ubt(I,j))
5232 elseif ((dt*ubt(i,j) < btcl_u(i,j)%uBT_EE) .and. (dt*uhbt(i,j) < btcl_u(i,j)%uh_EE)) then
5233 ! Expand the cubic fit to use this new point. ubt is negative.
5234 btcl_u(i,j)%ubt_EE = dt * ubt(i,j)
5235 if (3.0*uhbt(i,j) < 2.0*ubt(i,j) * btcl_u(i,j)%FA_u_E0) then
5236 ! No further bounding is needed.
5237 btcl_u(i,j)%uh_crvE = (uhbt(i,j) - ubt(i,j) * btcl_u(i,j)%FA_u_E0) / (dt**2 * ubt(i,j)**3)
5238 else ! This should not happen often!
5239 btcl_u(i,j)%FA_u_E0 = 1.5*uhbt(i,j) / ubt(i,j)
5240 btcl_u(i,j)%uh_crvE = -0.5*uhbt(i,j) / (dt**2 * ubt(i,j)**3)
5241 endif
5242 btcl_u(i,j)%uh_EE = dt * uhbt(i,j)
5243 ! I don't know whether this is helpful.
5244! BTCL_u(I,j)%FA_u_EE = min(BTCL_u(I,j)%FA_u_EE, uhbt(I,j) / ubt(I,j))
5245 endif
5246 enddo ; enddo
5247 !$OMP parallel do default(shared)
5248 do j=js-hs-1,je+hs ; do i=is-hs,ie+hs
5249 if ((dt*vbt(i,j) > btcl_v(i,j)%vBT_SS) .and. (dt*vhbt(i,j) > btcl_v(i,j)%vh_SS)) then
5250 ! Expand the cubic fit to use this new point. vbt is negative.
5251 btcl_v(i,j)%vbt_SS = dt * vbt(i,j)
5252 if (3.0*vhbt(i,j) < 2.0*vbt(i,j) * btcl_v(i,j)%FA_v_S0) then
5253 ! No further bounding is needed.
5254 btcl_v(i,j)%vh_crvS = (vhbt(i,j) - vbt(i,j) * btcl_v(i,j)%FA_v_S0) / (dt**2 * vbt(i,j)**3)
5255 else ! This should not happen often!
5256 btcl_v(i,j)%FA_v_S0 = 1.5*vhbt(i,j) / (vbt(i,j))
5257 btcl_v(i,j)%vh_crvS = -0.5*vhbt(i,j) / (dt**2 * vbt(i,j)**3)
5258 endif
5259 btcl_v(i,j)%vh_SS = dt * vhbt(i,j)
5260 ! I don't know whether this is helpful.
5261! BTCL_v(i,J)%FA_v_SS = min(BTCL_v(i,J)%FA_v_SS, vhbt(i,J) / vbt(i,J))
5262 elseif ((dt*vbt(i,j) < btcl_v(i,j)%vBT_NN) .and. (dt*vhbt(i,j) < btcl_v(i,j)%vh_NN)) then
5263 ! Expand the cubic fit to use this new point. vbt is negative.
5264 btcl_v(i,j)%vbt_NN = dt * vbt(i,j)
5265 if (3.0*vhbt(i,j) < 2.0*vbt(i,j) * btcl_v(i,j)%FA_v_N0) then
5266 ! No further bounding is needed.
5267 btcl_v(i,j)%vh_crvN = (vhbt(i,j) - vbt(i,j) * btcl_v(i,j)%FA_v_N0) / (dt**2 * vbt(i,j)**3)
5268 else ! This should not happen often!
5269 btcl_v(i,j)%FA_v_N0 = 1.5*vhbt(i,j) / (vbt(i,j))
5270 btcl_v(i,j)%vh_crvN = -0.5*vhbt(i,j) / (dt**2 * vbt(i,j)**3)
5271 endif
5272 btcl_v(i,j)%vh_NN = dt * vhbt(i,j)
5273 ! I don't know whether this is helpful.
5274! BTCL_v(i,J)%FA_v_NN = min(BTCL_v(i,J)%FA_v_NN, vhbt(i,J) / vbt(i,J))
5275 endif
5276 enddo ; enddo
5277
5278end subroutine adjust_local_bt_cont_types
5279
5280!> This subroutine uses the BT_cont_type to find the maximum face
5281!! areas, which can then be used for finding wave speeds, etc.
5282subroutine bt_cont_to_face_areas(BT_cont, Datu, Datv, G, US, MS, halo)
5283 type(bt_cont_type), intent(inout) :: BT_cont !< The BT_cont_type input to the
5284 !! barotropic solver.
5285 type(memory_size_type), intent(in) :: MS !< A type that describes the memory
5286 !! sizes of the argument arrays.
5287 real, dimension(MS%isdw-1:MS%iedw,MS%jsdw:MS%jedw), &
5288 intent(out) :: Datu !< The effective zonal face area [H L ~> m2 or kg m-1].
5289 real, dimension(MS%isdw:MS%iedw,MS%jsdw-1:MS%jedw), &
5290 intent(out) :: Datv !< The effective meridional face area [H L ~> m2 or kg m-1].
5291 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
5292 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
5293 integer, optional, intent(in) :: halo !< The extra halo size to use here.
5294
5295 ! Local variables
5296 integer :: i, j, is, ie, js, je, hs
5297 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
5298 hs = 1 ; if (present(halo)) hs = max(halo,0)
5299
5300 do j=js-hs,je+hs ; do i=is-1-hs,ie+hs
5301 datu(i,j) = max(bt_cont%FA_u_EE(i,j), bt_cont%FA_u_E0(i,j), &
5302 bt_cont%FA_u_W0(i,j), bt_cont%FA_u_WW(i,j))
5303 enddo ; enddo
5304 do j=js-1-hs,je+hs ; do i=is-hs,ie+hs
5305 datv(i,j) = max(bt_cont%FA_v_NN(i,j), bt_cont%FA_v_N0(i,j), &
5306 bt_cont%FA_v_S0(i,j), bt_cont%FA_v_SS(i,j))
5307 enddo ; enddo
5308
5309end subroutine bt_cont_to_face_areas
5310
5311!> Swap the values of two real variables
5312subroutine swap(a,b)
5313 real, intent(inout) :: a !< The first variable to be swapped [arbitrary units]
5314 real, intent(inout) :: b !< The second variable to be swapped [arbitrary units]
5315 real :: tmp ! A temporary variable [arbitrary units]
5316 tmp = a ; a = b ; b = tmp
5317end subroutine swap
5318
5319!> This subroutine determines the open face areas of cells for calculating
5320!! the barotropic transport.
5321subroutine find_face_areas(Datu, Datv, G, GV, US, CS, MS, halo, eta, add_max)
5322 type(memory_size_type), intent(in) :: MS !< A type that describes the memory sizes of the argument arrays.
5323 real, dimension(MS%isdw-1:MS%iedw,MS%jsdw:MS%jedw), &
5324 intent(out) :: Datu !< The open zonal face area [H L ~> m2 or kg m-1].
5325 real, dimension(MS%isdw:MS%iedw,MS%jsdw-1:MS%jedw), &
5326 intent(out) :: Datv !< The open meridional face area [H L ~> m2 or kg m-1].
5327 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
5328 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
5329 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
5330 type(barotropic_cs), intent(in) :: CS !< Barotropic control structure
5331 integer, intent(in) :: halo !< The halo size to use, default = 1.
5332 real, dimension(MS%isdw:MS%iedw,MS%jsdw:MS%jedw), &
5333 optional, intent(in) :: eta !< The barotropic free surface height anomaly
5334 !! or column mass anomaly [H ~> m or kg m-2].
5335 real, optional, intent(in) :: add_max !< A value to add to the maximum depth (used
5336 !! to overestimate the external wave speed) [Z ~> m].
5337
5338 ! Local variables
5339 real :: H1, H2 ! Temporary total thicknesses [H ~> m or kg m-2].
5340 real :: Z_to_H ! A local conversion factor [H Z-1 ~> nondim or kg m-3]
5341 integer :: i, j, is, ie, js, je, hs
5342 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
5343 hs = max(halo,0)
5344
5345 !$OMP parallel default(shared) private(H1,H2,Z_to_H)
5346 if (present(eta)) then
5347 ! The use of harmonic mean thicknesses ensure positive definiteness.
5348 if (gv%Boussinesq) then
5349 !$OMP do
5350 do j=js-hs,je+hs ; do i=is-1-hs,ie+hs
5351 h1 = cs%bathyT(i,j)*gv%Z_to_H + eta(i,j) ; h2 = cs%bathyT(i+1,j)*gv%Z_to_H + eta(i+1,j)
5352 datu(i,j) = 0.0 ; if ((h1 > 0.0) .and. (h2 > 0.0)) &
5353 datu(i,j) = cs%dy_Cu(i,j) * (2.0 * h1 * h2) / (h1 + h2)
5354 ! Datu(I,j) = CS%dy_Cu(I,j) * 0.5 * (H1 + H2)
5355 enddo ; enddo
5356 !$OMP do
5357 do j=js-1-hs,je+hs ; do i=is-hs,ie+hs
5358 h1 = cs%bathyT(i,j)*gv%Z_to_H + eta(i,j) ; h2 = cs%bathyT(i,j+1)*gv%Z_to_H + eta(i,j+1)
5359 datv(i,j) = 0.0 ; if ((h1 > 0.0) .and. (h2 > 0.0)) &
5360 datv(i,j) = cs%dx_Cv(i,j) * (2.0 * h1 * h2) / (h1 + h2)
5361 ! Datv(i,J) = CS%dy_v(i,J) * 0.5 * (H1 + H2)
5362 enddo ; enddo
5363 else
5364 !$OMP do
5365 do j=js-hs,je+hs ; do i=is-1-hs,ie+hs
5366 datu(i,j) = 0.0 ; if ((eta(i,j) > 0.0) .and. (eta(i+1,j) > 0.0)) &
5367 datu(i,j) = cs%dy_Cu(i,j) * (2.0 * eta(i,j) * eta(i+1,j)) / &
5368 (eta(i,j) + eta(i+1,j))
5369 ! Datu(I,j) = CS%dy_Cu(I,j) * 0.5 * (eta(i,j) + eta(i+1,j))
5370 enddo ; enddo
5371 !$OMP do
5372 do j=js-1-hs,je+hs ; do i=is-hs,ie+hs
5373 datv(i,j) = 0.0 ; if ((eta(i,j) > 0.0) .and. (eta(i,j+1) > 0.0)) &
5374 datv(i,j) = cs%dx_Cv(i,j) * (2.0 * eta(i,j) * eta(i,j+1)) / &
5375 (eta(i,j) + eta(i,j+1))
5376 ! Datv(i,J) = CS%dy_v(i,J) * 0.5 * (eta(i,j) + eta(i,j+1))
5377 enddo ; enddo
5378 endif
5379 elseif (present(add_max)) then
5380 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
5381
5382 !$OMP do
5383 do j=js-hs,je+hs ; do i=is-1-hs,ie+hs
5384 h1 = max((g%meanSL(i+1,j) + add_max) + g%bathyT(i+1,j), 0.0)
5385 h2 = max((g%meanSL(i,j) + add_max) + g%bathyT(i,j), 0.0)
5386 datu(i,j) = cs%dy_Cu(i,j) * z_to_h * max(h1, h2)
5387 enddo ; enddo
5388 !$OMP do
5389 do j=js-1-hs,je+hs ; do i=is-hs,ie+hs
5390 h1 = max((g%meanSL(i,j+1) + add_max) + g%bathyT(i,j+1), 0.0)
5391 h2 = max((g%meanSL(i,j) + add_max) + g%bathyT(i,j), 0.0)
5392 datv(i,j) = cs%dx_Cv(i,j) * z_to_h * max(h1, h2)
5393 enddo ; enddo
5394 else
5395 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
5396
5397 !$OMP do
5398 do j=js-hs,je+hs ; do i=is-1-hs,ie+hs
5399 h1 = max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) * z_to_h
5400 h2 = max(g%meanSL(i+1,j) + g%bathyT(i+1,j), 0.0) * z_to_h
5401 datu(i,j) = 0.0
5402 if ((h1 > 0.0) .and. (h2 > 0.0)) &
5403 datu(i,j) = cs%dy_Cu(i,j) * (2.0 * h1 * h2) / (h1 + h2)
5404 enddo ; enddo
5405 !$OMP do
5406 do j=js-1-hs,je+hs ; do i=is-hs,ie+hs
5407 h1 = max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) * z_to_h
5408 h2 = max(g%meanSL(i,j+1) + g%bathyT(i,j+1), 0.0) * z_to_h
5409 datv(i,j) = 0.0
5410 if ((h1 > 0.0) .and. (h2 > 0.0)) &
5411 datv(i,j) = cs%dx_Cv(i,j) * (2.0 * h1 * h2) / (h1 + h2)
5412 enddo ; enddo
5413 endif
5414 !$OMP end parallel
5415
5416end subroutine find_face_areas
5417
5418!> bt_mass_source determines the appropriately limited mass source for
5419!! the barotropic solver, along with a corrective fictitious mass source that
5420!! will drive the barotropic estimate of the free surface height toward the
5421!! baroclinic estimate.
5422subroutine bt_mass_source(h, eta, set_cor, G, GV, CS)
5423 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
5424 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
5425 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
5426 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: eta !< The free surface height that is to be
5427 !! corrected [H ~> m or kg m-2].
5428 logical, intent(in) :: set_cor !< A flag to indicate whether to set the corrective
5429 !! fluxes (and update the slowly varying part of eta_cor)
5430 !! (.true.) or whether to incrementally update the
5431 !! corrective fluxes.
5432 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
5433
5434 ! Local variables
5435 real :: h_tot(szi_(g)) ! The sum of the layer thicknesses [H ~> m or kg m-2].
5436 real :: eta_h(szi_(g)) ! The free surface height determined from
5437 ! the sum of the layer thicknesses [H ~> m or kg m-2].
5438 real :: d_eta ! The difference between estimates of the total
5439 ! thicknesses [H ~> m or kg m-2].
5440 integer :: is, ie, js, je, nz, i, j, k
5441
5442 if (.not.cs%module_is_initialized) call mom_error(fatal, "bt_mass_source: "// &
5443 "Module MOM_barotropic must be initialized before it is used.")
5444
5445 if (.not.cs%split) return
5446
5447 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
5448
5449 !$OMP parallel do default(shared) private(eta_h,h_tot,d_eta)
5450 do j=js,je
5451 do i=is,ie ; h_tot(i) = h(i,j,1) ; enddo
5452 if (gv%Boussinesq) then
5453 do i=is,ie ; eta_h(i) = h(i,j,1) - g%bathyT(i,j)*gv%Z_to_H ; enddo
5454 else
5455 do i=is,ie ; eta_h(i) = h(i,j,1) ; enddo
5456 endif
5457 do k=2,nz ; do i=is,ie
5458 eta_h(i) = eta_h(i) + h(i,j,k)
5459 h_tot(i) = h_tot(i) + h(i,j,k)
5460 enddo ; enddo
5461
5462 if (set_cor) then
5463 do i=is,ie
5464 d_eta = eta_h(i) - eta(i,j)
5465 cs%eta_cor(i,j) = d_eta
5466 enddo
5467 else
5468 do i=is,ie
5469 d_eta = eta_h(i) - eta(i,j)
5470 cs%eta_cor(i,j) = cs%eta_cor(i,j) + d_eta
5471 enddo
5472 endif
5473 enddo
5474
5475end subroutine bt_mass_source
5476
5477!> barotropic_init initializes a number of time-invariant fields used in the
5478!! barotropic calculation and initializes any barotropic fields that have not
5479!! already been initialized.
5480subroutine barotropic_init(u, v, h, Time, G, GV, US, param_file, diag, CS, &
5481 restart_CS, calc_dtbt, BT_cont, OBC, SAL_CSp, HA_CSp)
5482 type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
5483 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
5484 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
5485 real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
5486 intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
5487 real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
5488 intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
5489 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
5490 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
5491 type(time_type), target, intent(in) :: time !< The current model time.
5492 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
5493 type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic
5494 !! output.
5495 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
5496 type(mom_restart_cs), intent(in) :: restart_cs !< MOM restart control structure
5497 logical, intent(out) :: calc_dtbt !< If true, the barotropic time step must
5498 !! be recalculated before stepping.
5499 type(bt_cont_type), pointer :: bt_cont !< A structure with elements that describe the
5500 !! effective open face areas as a function of
5501 !! barotropic flow.
5502 type(ocean_obc_type), pointer :: obc !< The open boundary condition structure.
5503 type(sal_cs), target, optional :: sal_csp !< A pointer to the control structure of the
5504 !! SAL module.
5505 type(harmonic_analysis_cs), target, optional :: ha_csp !< A pointer to the control structure of the
5506 !! harmonic analysis module
5507
5508 ! This include declares and sets the variable "version".
5509# include "version_variable.h"
5510 ! Local variables
5511 character(len=40) :: mdl = "MOM_barotropic" ! This module's name.
5512 real :: datu(szibs_(g),szj_(g)) ! Zonal open face area [H L ~> m2 or kg m-1].
5513 real :: datv(szi_(g),szjbs_(g)) ! Meridional open face area [H L ~> m2 or kg m-1].
5514 real :: gtot_estimate ! Summed GV%g_prime [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2], to give an
5515 ! upper-bound estimate for pbce.
5516 real :: ssh_extra ! An estimate of how much higher SSH might get, for use
5517 ! in calculating the safe external wave speed [Z ~> m].
5518 real :: dtbt_input ! The input value of DTBT, [nondim] if negative or [s] if positive.
5519 real :: dtbt_restart ! A temporary copy of CS%dtbt read from a restart file [T ~> s]
5520 real :: wave_drag_scale ! A scaling factor for the barotropic linear wave drag
5521 ! piston velocities [nondim].
5522 character(len=200) :: inputdir ! The directory in which to find input files.
5523 character(len=200) :: wave_drag_file ! The file from which to read the wave
5524 ! drag piston velocity.
5525 character(len=80) :: wave_drag_var ! The wave drag piston velocity variable
5526 ! name in wave_drag_file.
5527 character(len=80) :: wave_drag_u ! The wave drag piston velocity variable
5528 ! name in wave_drag_file.
5529 character(len=80) :: wave_drag_v ! The wave drag piston velocity variable
5530 ! name in wave_drag_file.
5531 real :: htot ! Total column thickness used when BT_NONLIN_STRESS is false [Z ~> m].
5532 real :: z_to_h ! A local unit conversion factor [H Z-1 ~> nondim or kg m-3]
5533 real :: h_to_z ! A local unit conversion factor [Z H-1 ~> nondim or m3 kg-1]
5534 real :: det_de ! The partial derivative due to self-attraction and loading of the reference
5535 ! geopotential with the sea surface height when scalar SAL are enabled [nondim].
5536 ! This is typically ~0.09 or less.
5537 real :: h_a_neglect ! A cell volume or mass that is so small it is usually lost
5538 ! in roundoff and can be neglected [H L2 ~> m3 or kg]
5539 real, allocatable :: lin_drag_h(:,:) ! A spatially varying linear drag coefficient at tracer points
5540 ! that acts on the barotropic flow [H T-1 ~> m s-1 or kg m-2 s-1].
5541
5542 type(memory_size_type) :: ms
5543 type(group_pass_type) :: pass_static_data, pass_q_d_cor
5544 type(group_pass_type) :: pass_bt_hbt_btav, pass_a_polarity
5545 integer :: default_answer_date ! The default setting for the various ANSWER_DATE flags.
5546 logical :: use_bt_cont_type
5547 logical :: mask_coastal_pressure_force ! If true, apply masks to some stored inverse grid spacings
5548 ! so that diagnosed barotropic pressure gradient forces are zero at
5549 ! land, coastal or OBC points.
5550 logical :: use_tides
5551 logical :: obc_projection_bug
5552 logical :: enable_bugs ! If true, the defaults for recently added bug-fix flags are set to
5553 ! recreate the bugs, or if false bugs are only used if actively selected.
5554 logical :: visc_rem_bug ! Stores the value of runtime paramter VISC_REM_BUG.
5555 character(len=48) :: thickness_units, flux_units
5556 character*(40) :: hvel_str
5557 integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
5558 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
5559 integer :: isdw, iedw, jsdw, jedw
5560 integer :: i, j, k
5561 integer :: wd_halos(2), bt_halo_sz
5562 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
5563 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
5564 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
5565 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
5566 ms%isdw = g%isd ; ms%iedw = g%ied ; ms%jsdw = g%jsd ; ms%jedw = g%jed
5567
5568 if (cs%module_is_initialized) then
5569 call mom_error(warning, "barotropic_init called with a control structure "// &
5570 "that has already been initialized.")
5571 return
5572 endif
5573 cs%module_is_initialized = .true.
5574
5575 cs%diag => diag ; cs%Time => time
5576 if (present(sal_csp)) then
5577 cs%SAL_CSp => sal_csp
5578 endif
5579
5580 ! Read all relevant parameters and write them to the model log.
5581 call get_param(param_file, mdl, "SPLIT", cs%split, default=.true., do_not_log=.true.)
5582 call log_version(param_file, mdl, version, "", log_to_all=.true., layout=cs%split, &
5583 debugging=cs%split, all_default=.not.cs%split)
5584 call get_param(param_file, mdl, "SPLIT", cs%split, &
5585 "Use the split time stepping if true.", default=.true.)
5586 if (.not.cs%split) return
5587
5588 call get_param(param_file, mdl, "USE_BT_CONT_TYPE", use_bt_cont_type, &
5589 "If true, use a structure with elements that describe "//&
5590 "effective face areas from the summed continuity solver "//&
5591 "as a function the barotropic flow in coupling between "//&
5592 "the barotropic and baroclinic flow. This is only used "//&
5593 "if SPLIT is true.", default=.true.)
5594 call get_param(param_file, mdl, "INTEGRAL_BT_CONTINUITY", cs%integral_bt_cont, &
5595 "If true, use the time-integrated velocity over the barotropic steps "//&
5596 "to determine the integrated transports used to update the continuity "//&
5597 "equation. Otherwise the transports are the sum of the transports based on "//&
5598 "a series of instantaneous velocities and the BT_CONT_TYPE for transports. "//&
5599 "This is only valid if USE_BT_CONT_TYPE = True.", &
5600 default=.false., do_not_log=.not.use_bt_cont_type)
5601 call get_param(param_file, mdl, "BOUND_BT_CORRECTION", cs%bound_BT_corr, &
5602 "If true, the corrective pseudo mass-fluxes into the "//&
5603 "barotropic solver are limited to values that require "//&
5604 "less than maxCFL_BT_cont to be accommodated.",default=.false.)
5605 call get_param(param_file, mdl, "BT_CONT_CORR_BOUNDS", cs%BT_cont_bounds, &
5606 "If true, and BOUND_BT_CORRECTION is true, use the "//&
5607 "BT_cont_type variables to set limits determined by "//&
5608 "MAXCFL_BT_CONT on the CFL number of the velocities "//&
5609 "that are likely to be driven by the corrective mass fluxes.", &
5610 default=.true., do_not_log=.not.cs%bound_BT_corr)
5611 call get_param(param_file, mdl, "ADJUST_BT_CONT", cs%adjust_BT_cont, &
5612 "If true, adjust the curve fit to the BT_cont type "//&
5613 "that is used by the barotropic solver to match the "//&
5614 "transport about which the flow is being linearized.", &
5615 default=.false., do_not_log=.not.use_bt_cont_type)
5616 call get_param(param_file, mdl, "GRADUAL_BT_ICS", cs%gradual_BT_ICs, &
5617 "If true, adjust the initial conditions for the "//&
5618 "barotropic solver to the values from the layered "//&
5619 "solution over a whole timestep instead of instantly. "//&
5620 "This is a decent approximation to the inclusion of "//&
5621 "sum(u dh_dt) while also correcting for truncation errors.", &
5622 default=.false.)
5623 call get_param(param_file, mdl, "BT_ADJUST_SRC_FOR_FILTER", cs%bt_adjust_src_for_filter, &
5624 "If true, increases the rate at which BT mass sources are applied so "//&
5625 "that they are all used up before the filtering period starts. "//&
5626 "This option is only valid if INTEGRAL_BT_CONTINUITY = True.", &
5627 default=.false., do_not_log=.not.cs%integral_bt_cont)
5628 call get_param(param_file, mdl, "BT_LIMIT_INTEGRAL_TRANSPORT", cs%bt_limit_integral_transport, &
5629 "If true, limit the time-integrated transports by the initial volume "//&
5630 "accounting for sinks of mass. The limiter uses MAXCFL_BT_CONT. "//&
5631 "This option is only valid if INTEGRAL_BT_CONTINUITY = True.", &
5632 default=.false., do_not_log=.not.cs%integral_bt_cont)
5633 call get_param(param_file, mdl, "BT_USE_VISC_REM_U_UH0", cs%visc_rem_u_uh0, &
5634 "If true, use the viscous remnants when estimating the "//&
5635 "barotropic velocities that were used to calculate uh0 "//&
5636 "and vh0. False is probably the better choice.", default=.false.)
5637 call get_param(param_file, mdl, "BT_USE_WIDE_HALOS", cs%use_wide_halos, &
5638 "If true, use wide halos and march in during the "//&
5639 "barotropic time stepping for efficiency.", default=.true., &
5640 layoutparam=.true.)
5641 call get_param(param_file, mdl, "BTHALO", bt_halo_sz, &
5642 "The minimum halo size for the barotropic solver.", default=0, &
5643 layoutparam=.true.)
5644 call get_param(param_file, mdl, "BT_WIDE_HALO_MIN_STENCIL", cs%min_stencil, &
5645 "The minimum stencil width to use with the wide halo iterations. "//&
5646 "A nonzero value may be useful for debugging purposes, but at the "//&
5647 "cost of reducing the efficiency gain from BT_USE_WIDE_HALOS.", &
5648 default=0, layoutparam=.true., do_not_log=.not.cs%use_wide_halos)
5649#ifdef STATIC_MEMORY_
5650 if ((bt_halo_sz > 0) .and. (bt_halo_sz /= bthalo_)) call mom_error(fatal, &
5651 "barotropic_init: Run-time values of BTHALO must agree with the "//&
5652 "macro BTHALO_ with STATIC_MEMORY_.")
5653 wd_halos(1) = whaloi_+nihalo_ ; wd_halos(2) = whaloj_+njhalo_
5654#else
5655 wd_halos(1) = bt_halo_sz ; wd_halos(2) = bt_halo_sz
5656#endif
5657 call get_param(param_file, mdl, "NONLINEAR_BT_CONTINUITY", cs%Nonlinear_continuity, &
5658 "If true, use nonlinear transports in the barotropic "//&
5659 "continuity equation. This does not apply if "//&
5660 "USE_BT_CONT_TYPE is true.", default=.false., do_not_log=use_bt_cont_type)
5661 call get_param(param_file, mdl, "NONLIN_BT_CONT_UPDATE_PERIOD", cs%Nonlin_cont_update_period, &
5662 "If NONLINEAR_BT_CONTINUITY is true, this is the number "//&
5663 "of barotropic time steps between updates to the face "//&
5664 "areas, or 0 to update only before the barotropic stepping.", &
5665 default=1, do_not_log=.not.cs%Nonlinear_continuity)
5666
5667 call get_param(param_file, mdl, "BT_PROJECT_VELOCITY", cs%BT_project_velocity,&
5668 "If true, step the barotropic velocity first and project "//&
5669 "out the velocity tendency by 1+BEBT when calculating the "//&
5670 "transport. The default (false) is to use a predictor "//&
5671 "continuity step to find the pressure field, and then "//&
5672 "to do a corrector continuity step using a weighted "//&
5673 "average of the old and new velocities, with weights "//&
5674 "of (1-BEBT) and BEBT.", default=.false.)
5675 call get_param(param_file, mdl, "BT_NONLIN_STRESS", cs%nonlin_stress, &
5676 "If true, use the full depth of the ocean at the start of the barotropic "//&
5677 "step when calculating the surface stress contribution to the barotropic "//&
5678 "acclerations. Otherwise use the depth based on bathyT.", default=.false.)
5679 call get_param(param_file, mdl, "BT_RHO_LINEARIZED", cs%Rho_BT_lin, &
5680 "A density that is used to convert total water column thicknesses into mass "//&
5681 "in non-Boussinesq mode with linearized options in the barotropic solver or "//&
5682 "when estimating the stable barotropic timestep without access to the full "//&
5683 "baroclinic model state.", &
5684 units="kg m-3", default=gv%Rho0*us%R_to_kg_m3, scale=us%kg_m3_to_R, &
5685 do_not_log=gv%Boussinesq)
5686
5687 call get_param(param_file, mdl, "DYNAMIC_SURFACE_PRESSURE", cs%dynamic_psurf, &
5688 "If true, add a dynamic pressure due to a viscous ice "//&
5689 "shelf, for instance.", default=.false.)
5690 call get_param(param_file, mdl, "ICE_LENGTH_DYN_PSURF", cs%ice_strength_length, &
5691 "The length scale at which the Rayleigh damping rate due "//&
5692 "to the ice strength should be the same as if a Laplacian "//&
5693 "were applied, if DYNAMIC_SURFACE_PRESSURE is true.", &
5694 units="m", default=1.0e4, scale=us%m_to_L, do_not_log=.not.cs%dynamic_psurf)
5695 call get_param(param_file, mdl, "DEPTH_MIN_DYN_PSURF", cs%Dmin_dyn_psurf, &
5696 "The minimum depth to use in limiting the size of the "//&
5697 "dynamic surface pressure for stability, if "//&
5698 "DYNAMIC_SURFACE_PRESSURE is true..", &
5699 units="m", default=1.0e-6, scale=gv%m_to_H, do_not_log=.not.cs%dynamic_psurf)
5700 call get_param(param_file, mdl, "CONST_DYN_PSURF", cs%const_dyn_psurf, &
5701 "The constant that scales the dynamic surface pressure, "//&
5702 "if DYNAMIC_SURFACE_PRESSURE is true. Stable values "//&
5703 "are < ~1.0.", units="nondim", default=0.9, do_not_log=.not.cs%dynamic_psurf)
5704
5705 call get_param(param_file, mdl, "BT_CORIOLIS_SCALE", cs%BT_Coriolis_scale, &
5706 "A factor by which the barotropic Coriolis anomaly terms are scaled.", &
5707 units="nondim", default=1.0)
5708 call get_param(param_file, mdl, "DEFAULT_ANSWER_DATE", default_answer_date, &
5709 "This sets the default value for the various _ANSWER_DATE parameters.", &
5710 default=99991231)
5711 call get_param(param_file, mdl, "BAROTROPIC_ANSWER_DATE", cs%answer_date, &
5712 "The vintage of the expressions in the barotropic solver. "//&
5713 "Values below 20190101 recover the answers from the end of 2018, "//&
5714 "while higher values use more efficient or general expressions.", &
5715 default=default_answer_date, do_not_log=.not.gv%Boussinesq)
5716 if (.not.gv%Boussinesq) cs%answer_date = max(cs%answer_date, 20230701)
5717
5718 call get_param(param_file, mdl, "ENABLE_BUGS_BY_DEFAULT", enable_bugs, &
5719 default=.true., do_not_log=.true.) ! This is logged from MOM.F90.
5720 call get_param(param_file, mdl, "VISC_REM_BUG", visc_rem_bug, default=.false., do_not_log=.true.)
5721 call get_param(param_file, mdl, "VISC_REM_BT_WEIGHT_BUG", cs%wt_uv_bug, &
5722 "If true, recover a bug in barotropic solver that uses an unnormalized weight "//&
5723 "function for vertical averages of baroclinic velocity and forcing. Default "//&
5724 "of this flag is set by VISC_REM_BUG.", default=visc_rem_bug)
5725 call get_param(param_file, mdl, "EXTERIOR_OBC_BUG", cs%exterior_OBC_bug, &
5726 "If true, recover a bug in barotropic solver and other routines when "//&
5727 "boundary contitions interior to the domain are used.", &
5728 default=enable_bugs, do_not_log=.true.)
5729 call get_param(param_file, mdl, "OBC_PROJECTION_BUG", obc_projection_bug, &
5730 "If false, use only interior ocean points at OBCs to specify several "//&
5731 "calculations at OBC points, and it avoids applying a land mask at the bay-like "//&
5732 "intersection of orthogonal OBC segments. Otherwise the calculation of terms "//&
5733 "like the potential vorticity used in the barotropic solver relies on bathymetry "//&
5734 "or other fields being projected outward across OBCs. This option changes "//&
5735 "answers for some configurations that use OBCs.", &
5736 default=enable_bugs, do_not_log=.true.)
5737 cs%interior_OBC_PV = .not.obc_projection_bug
5738
5739 call get_param(param_file, mdl, "TIDES", use_tides, &
5740 "If true, apply tidal momentum forcing.", default=.false.)
5741 if (use_tides .and. present(ha_csp)) cs%HA_CSp => ha_csp
5742 call get_param(param_file, mdl, "CALCULATE_SAL", cs%calculate_SAL, &
5743 "If true, calculate self-attraction and loading.", default=use_tides)
5744 det_de = 0.0
5745 if (cs%calculate_SAL .and. associated(cs%SAL_CSp)) &
5746 call scalar_sal_sensitivity(cs%SAL_CSp, det_de)
5747 call get_param(param_file, mdl, "BAROTROPIC_TIDAL_SAL_BUG", cs%tidal_sal_bug, &
5748 "If true, the tidal self-attraction and loading anomaly in the barotropic "//&
5749 "solver has the wrong sign, replicating a long-standing bug with a scalar "//&
5750 "self-attraction and loading term or the SAL term from a previous simulation.", &
5751 default=.false., do_not_log=(det_de==0.0))
5752 call get_param(param_file, mdl, "TIDAL_SAL_FLATHER", cs%tidal_sal_flather, &
5753 "If true, then apply adjustments to the external gravity "//&
5754 "wave speed used with the Flather OBC routine consistent "//&
5755 "with the barotropic solver. This applies to cases with "//&
5756 "tidal forcing using the scalar self-attraction approximation. "//&
5757 "The default is currently False in order to retain previous answers "//&
5758 "but should be set to True for new experiments", default=.false.)
5759
5760 call get_param(param_file, mdl, "SADOURNY", cs%Sadourny, &
5761 "If true, the Coriolis terms are discretized with the "//&
5762 "Sadourny (1975) energy conserving scheme, otherwise "//&
5763 "the Arakawa & Hsu scheme is used. If the internal "//&
5764 "deformation radius is not resolved, the Sadourny scheme "//&
5765 "should probably be used.", default=.true.)
5766
5767 call get_param(param_file, mdl, "BT_THICK_SCHEME", hvel_str, &
5768 "A string describing the scheme that is used to set the "//&
5769 "open face areas used for barotropic transport and the "//&
5770 "relative weights of the accelerations. Valid values are:\n"//&
5771 "\t ARITHMETIC - arithmetic mean layer thicknesses \n"//&
5772 "\t HARMONIC - harmonic mean layer thicknesses \n"//&
5773 "\t HYBRID (the default) - use arithmetic means for \n"//&
5774 "\t layers above the shallowest bottom, the harmonic \n"//&
5775 "\t mean for layers below, and a weighted average for \n"//&
5776 "\t layers that straddle that depth \n"//&
5777 "\t FROM_BT_CONT - use the average thicknesses kept \n"//&
5778 "\t in the h_u and h_v fields of the BT_cont_type", &
5779 default=bt_cont_string)
5780 select case (hvel_str)
5781 case (hybrid_string) ; cs%hvel_scheme = hybrid
5782 case (harmonic_string) ; cs%hvel_scheme = harmonic
5783 case (arithmetic_string) ; cs%hvel_scheme = arithmetic
5784 case (bt_cont_string) ; cs%hvel_scheme = from_bt_cont
5785 case default
5786 call mom_mesg('barotropic_init: BT_THICK_SCHEME ="'//trim(hvel_str)//'"', 0)
5787 call mom_error(fatal, "barotropic_init: Unrecognized setting "// &
5788 "#define BT_THICK_SCHEME "//trim(hvel_str)//" found in input file.")
5789 end select
5790 if ((cs%hvel_scheme == from_bt_cont) .and. .not.use_bt_cont_type) &
5791 call mom_error(fatal, "barotropic_init: BT_THICK_SCHEME FROM_BT_CONT "//&
5792 "can only be used if USE_BT_CONT_TYPE is defined.")
5793
5794 call get_param(param_file, mdl, "BT_STRONG_DRAG", cs%strong_drag, &
5795 "If true, use a stronger estimate of the retarding "//&
5796 "effects of strong bottom drag, by making it implicit "//&
5797 "with the barotropic time-step instead of implicit with "//&
5798 "the baroclinic time-step and dividing by the number of "//&
5799 "barotropic steps.", default=.false.)
5800 call get_param(param_file, mdl, "RESCALE_STRONG_DRAG", cs%rescale_strong_drag, &
5801 "If true, reduce the barotropic contribution to the layer accelerations "//&
5802 "to account for the difference between the forces that can be counteracted "//&
5803 "by the stronger drag with BT_STRONG_DRAG and the average of the layer "//&
5804 "viscous remnants after a baroclinic timestep.", &
5805 default=.false., do_not_log=.not.cs%strong_drag)
5806 call get_param(param_file, mdl, "BT_LINEAR_WAVE_DRAG", cs%linear_wave_drag, &
5807 "If true, apply a linear drag to the barotropic velocities, "//&
5808 "using rates set by lin_drag_u & _v divided by the depth of "//&
5809 "the ocean. This was introduced to facilitate tide modeling.", &
5810 default=.false.)
5811 call get_param(param_file, mdl, "BT_LINEAR_FREQ_DRAG", cs%linear_freq_drag, &
5812 "If true, apply frequency-dependent drag to the tidal velocities. "//&
5813 "The streaming band-pass filter must be turned on.", default=.false.)
5814 call get_param(param_file, mdl, "BT_WAVE_DRAG_FILE", wave_drag_file, &
5815 "The name of the file with the barotropic linear wave drag "//&
5816 "piston velocities.", default="", &
5817 do_not_log=(.not.cs%linear_wave_drag) .and. (.not.cs%linear_freq_drag))
5818 call get_param(param_file, mdl, "BT_WAVE_DRAG_VAR", wave_drag_var, &
5819 "The name of the variable in BT_WAVE_DRAG_FILE with the "//&
5820 "barotropic linear wave drag piston velocities at h points. "//&
5821 "It will not be used if both BT_WAVE_DRAG_U and BT_WAVE_DRAG_V are defined.", &
5822 default="rH", do_not_log=.not.cs%linear_wave_drag)
5823 call get_param(param_file, mdl, "BT_WAVE_DRAG_U", wave_drag_u, &
5824 "The name of the variable in BT_WAVE_DRAG_FILE with the "//&
5825 "barotropic linear wave drag piston velocities at u points.", &
5826 default="", do_not_log=.not.cs%linear_wave_drag)
5827 call get_param(param_file, mdl, "BT_WAVE_DRAG_V", wave_drag_v, &
5828 "The name of the variable in BT_WAVE_DRAG_FILE with the "//&
5829 "barotropic linear wave drag piston velocities at v points.", &
5830 default="", do_not_log=.not.cs%linear_wave_drag)
5831 call get_param(param_file, mdl, "BT_WAVE_DRAG_SCALE", wave_drag_scale, &
5832 "A scaling factor for the barotropic linear wave drag "//&
5833 "piston velocities.", default=1.0, units="nondim", &
5834 do_not_log=.not.cs%linear_wave_drag)
5835
5836 call get_param(param_file, mdl, "CLIP_BT_VELOCITY", cs%clip_velocity, &
5837 "If true, limit any velocity components that exceed "//&
5838 "CFL_TRUNCATE. This should only be used as a desperate "//&
5839 "debugging measure.", default=.false.)
5840 call get_param(param_file, mdl, "CFL_TRUNCATE", cs%CFL_trunc, &
5841 "The value of the CFL number that will cause velocity "//&
5842 "components to be truncated; instability can occur past 0.5.", &
5843 units="nondim", default=0.5, do_not_log=.not.cs%clip_velocity)
5844 call get_param(param_file, mdl, "MAXVEL", cs%maxvel, &
5845 "The maximum velocity allowed before the velocity "//&
5846 "components are truncated.", units="m s-1", default=3.0e8, scale=us%m_s_to_L_T, &
5847 do_not_log=.not.cs%clip_velocity)
5848 call get_param(param_file, mdl, "MAXCFL_BT_CONT", cs%maxCFL_BT_cont, &
5849 "The maximum permitted CFL number associated with the "//&
5850 "barotropic accelerations from the summed velocities "//&
5851 "times the time-derivatives of thicknesses.", units="nondim", &
5852 default=0.25)
5853 call get_param(param_file, mdl, "VEL_UNDERFLOW", cs%vel_underflow, &
5854 "A negligibly small velocity magnitude below which velocity "//&
5855 "components are set to 0. A reasonable value might be "//&
5856 "1e-30 m/s, which is less than an Angstrom divided by "//&
5857 "the age of the universe.", units="m s-1", default=0.0, scale=us%m_s_to_L_T)
5858
5859 call get_param(param_file, mdl, "DT_BT_FILTER", cs%dt_bt_filter, &
5860 "A time-scale over which the barotropic mode solutions "//&
5861 "are filtered, in seconds if positive, or as a fraction "//&
5862 "of DT if negative. When used this can never be taken to "//&
5863 "be longer than 2*dt. Set this to 0 to apply no filtering.", &
5864 units="sec or nondim", default=-0.25)
5865 if (cs%dt_bt_filter > 0.0) cs%dt_bt_filter = us%s_to_T*cs%dt_bt_filter
5866 call get_param(param_file, mdl, "G_BT_EXTRA", cs%G_extra, &
5867 "A nondimensional factor by which gtot is enhanced.", &
5868 units="nondim", default=0.0)
5869 call get_param(param_file, mdl, "SSH_EXTRA", ssh_extra, &
5870 "An estimate of how much higher SSH might get, for use "//&
5871 "in calculating the safe external wave speed. The "//&
5872 "default is the minimum of 10 m or 5% of MAXIMUM_DEPTH.", &
5873 units="m", default=min(10.0,0.05*g%max_depth*us%Z_to_m), scale=us%m_to_Z)
5874
5875 call get_param(param_file, mdl, "DEBUG", cs%debug, &
5876 "If true, write out verbose debugging data.", &
5877 default=.false., debuggingparam=.true.)
5878 call get_param(param_file, mdl, "DEBUG_BT", cs%debug_bt, &
5879 "If true, write out verbose debugging data within the "//&
5880 "barotropic time-stepping loop. The data volume can be "//&
5881 "quite large if this is true.", default=cs%debug, &
5882 debuggingparam=.true.)
5883 call get_param(param_file, mdl, "DEBUG_BT_WIDE_HALOS", cs%debug_wide_halos, &
5884 "If true, write the checksums on the full wide halos. Otherwise only the "//&
5885 "output for the final computational domain is written. This can be valuable "//&
5886 "for debugging certain cases where the stencil used in the wide halo "//&
5887 "iterations depends on which opoen boundary conditions are in the halos.", &
5888 default=.true., do_not_log=.not.(cs%debug_bt.and.cs%use_wide_halos), debuggingparam=.true.)
5889
5890 call get_param(param_file, mdl, "LINEARIZED_BT_CORIOLIS", cs%linearized_BT_PV, &
5891 "If true use the bottom depth instead of the total water column thickness "//&
5892 "in the barotropic Coriolis term calculations.", default=.true.)
5893 call get_param(param_file, mdl, "BEBT", cs%bebt, &
5894 "BEBT determines whether the barotropic time stepping "//&
5895 "uses the forward-backward time-stepping scheme or a "//&
5896 "backward Euler scheme. BEBT is valid in the range from "//&
5897 "0 (for a forward-backward treatment of nonrotating "//&
5898 "gravity waves) to 1 (for a backward Euler treatment). "//&
5899 "In practice, BEBT must be greater than about 0.05.", &
5900 units="nondim", default=0.1)
5901 ! Note that dtbt_input is not rescaled because it has different units for
5902 ! positive [s] and negative [nondim] values.
5903 call get_param(param_file, mdl, "DTBT", dtbt_input, &
5904 "The barotropic time step, in s. DTBT is only used with "//&
5905 "the split explicit time stepping. To set the time step "//&
5906 "automatically based the maximum stable value use 0, or "//&
5907 "a negative value gives the fraction of the stable value. "//&
5908 "Setting DTBT to 0 is the same as setting it to -0.98. "//&
5909 "The value of DTBT that will actually be used is an "//&
5910 "integer fraction of DT, rounding down.", &
5911 units="s or nondim", default=-0.98)
5912 call get_param(param_file, mdl, "BT_USE_OLD_CORIOLIS_BRACKET_BUG", cs%use_old_coriolis_bracket_bug, &
5913 "If True, use an order of operations that is not bitwise "//&
5914 "rotationally symmetric in the meridional Coriolis term of "//&
5915 "the barotropic solver.", default=.false.)
5916 call get_param(param_file, mdl, "MASK_COASTAL_PRESSURE_FORCE", mask_coastal_pressure_force, &
5917 "If true, use the land masks to zero out the diagnosed barotropic pressure "//&
5918 "gradient accelerations at coastal or land points. This changes diagnostics "//&
5919 "and improves the reproducibility of certain debugging checksums, but it "//&
5920 "does not alter the solutions themselves.", default=.false.)
5921 !### Change the default for MASK_COASTAL_PRESSURE_FORCE to true?
5922
5923 ! Initialize a version of the MOM domain that is specific to the barotropic solver.
5924 call clone_mom_domain(g%Domain, cs%BT_Domain, min_halo=wd_halos, symmetric=.true.)
5925#ifdef STATIC_MEMORY_
5926 if (wd_halos(1) /= whaloi_+nihalo_) call mom_error(fatal, "barotropic_init: "//&
5927 "Barotropic x-halo sizes are incorrectly resized with STATIC_MEMORY_.")
5928 if (wd_halos(2) /= whaloj_+njhalo_) call mom_error(fatal, "barotropic_init: "//&
5929 "Barotropic y-halo sizes are incorrectly resized with STATIC_MEMORY_.")
5930#else
5931 if (bt_halo_sz > 0) then
5932 if (wd_halos(1) > bt_halo_sz) &
5933 call mom_mesg("barotropic_init: barotropic x-halo size increased.", 3)
5934 if (wd_halos(2) > bt_halo_sz) &
5935 call mom_mesg("barotropic_init: barotropic y-halo size increased.", 3)
5936 endif
5937#endif
5938 call log_param(param_file, mdl, "!BT x-halo", wd_halos(1), &
5939 "The barotropic x-halo size that is actually used.", &
5940 layoutparam=.true.)
5941 call log_param(param_file, mdl, "!BT y-halo", wd_halos(2), &
5942 "The barotropic y-halo size that is actually used.", &
5943 layoutparam=.true.)
5944
5945 cs%isdw = g%isc-wd_halos(1) ; cs%iedw = g%iec+wd_halos(1)
5946 cs%jsdw = g%jsc-wd_halos(2) ; cs%jedw = g%jec+wd_halos(2)
5947 isdw = cs%isdw ; iedw = cs%iedw ; jsdw = cs%jsdw ; jedw = cs%jedw
5948
5949 alloc_(cs%frhatu(isdb:iedb,jsd:jed,nz)) ; alloc_(cs%frhatv(isd:ied,jsdb:jedb,nz))
5950 alloc_(cs%eta_cor(isd:ied,jsd:jed))
5951 if (cs%bound_BT_corr) &
5952 allocate(cs%eta_cor_bound(isd:ied,jsd:jed), source=0.0)
5953 alloc_(cs%IDatu(isdb:iedb,jsd:jed)) ; alloc_(cs%IDatv(isd:ied,jsdb:jedb))
5954
5955 alloc_(cs%ua_polarity(isdw:iedw,jsdw:jedw))
5956 alloc_(cs%va_polarity(isdw:iedw,jsdw:jedw))
5957
5958 cs%frhatu(:,:,:) = 0.0 ; cs%frhatv(:,:,:) = 0.0
5959 cs%eta_cor(:,:) = 0.0
5960 cs%IDatu(:,:) = 0.0 ; cs%IDatv(:,:) = 0.0
5961
5962 cs%ua_polarity(:,:) = 1.0 ; cs%va_polarity(:,:) = 1.0
5963 call create_group_pass(pass_a_polarity, cs%ua_polarity, cs%va_polarity, cs%BT_domain, to_all, agrid)
5964 call do_group_pass(pass_a_polarity, cs%BT_domain)
5965
5966 if (use_bt_cont_type) &
5967 call alloc_bt_cont_type(bt_cont, g, gv, (cs%hvel_scheme == from_bt_cont))
5968
5969 if (cs%debug) then ! Make a local copy of loop ranges for chksum calls
5970 allocate(cs%debug_BT_HI)
5971 cs%debug_BT_HI%isc = g%isc
5972 cs%debug_BT_HI%iec = g%iec
5973 cs%debug_BT_HI%jsc = g%jsc
5974 cs%debug_BT_HI%jec = g%jec
5975 cs%debug_BT_HI%IscB = g%isc-1
5976 cs%debug_BT_HI%IecB = g%iec
5977 cs%debug_BT_HI%JscB = g%jsc-1
5978 cs%debug_BT_HI%JecB = g%jec
5979 cs%debug_BT_HI%isd = cs%isdw
5980 cs%debug_BT_HI%ied = cs%iedw
5981 cs%debug_BT_HI%jsd = cs%jsdw
5982 cs%debug_BT_HI%jed = cs%jedw
5983 cs%debug_BT_HI%IsdB = cs%isdw-1
5984 cs%debug_BT_HI%IedB = cs%iedw
5985 cs%debug_BT_HI%JsdB = cs%jsdw-1
5986 cs%debug_BT_HI%JedB = cs%jedw
5987 cs%debug_BT_HI%turns = g%HI%turns
5988 endif
5989
5990 ! IareaT, IdxCu, and IdyCv need to be allocated with wide halos.
5991 alloc_(cs%IareaT(cs%isdw:cs%iedw,cs%jsdw:cs%jedw)) ; cs%IareaT(:,:) = 0.0
5992 alloc_(cs%bathyT(cs%isdw:cs%iedw,cs%jsdw:cs%jedw)) ; cs%bathyT(:,:) = 0.0
5993 alloc_(cs%IdxCu(cs%isdw-1:cs%iedw,cs%jsdw:cs%jedw)) ; cs%IdxCu(:,:) = 0.0
5994 alloc_(cs%IdyCv(cs%isdw:cs%iedw,cs%jsdw-1:cs%jedw)) ; cs%IdyCv(:,:) = 0.0
5995 alloc_(cs%dy_Cu(cs%isdw-1:cs%iedw,cs%jsdw:cs%jedw)) ; cs%dy_Cu(:,:) = 0.0
5996 alloc_(cs%dx_Cv(cs%isdw:cs%iedw,cs%jsdw-1:cs%jedw)) ; cs%dx_Cv(:,:) = 0.0
5997 allocate(cs%IareaT_OBCmask(isdw:iedw,jsdw:jedw), source=0.0)
5998 alloc_(cs%OBCmask_u(cs%isdw-1:cs%iedw,cs%jsdw:cs%jedw)) ; cs%OBCmask_u(:,:) = 0.0
5999 alloc_(cs%OBCmask_v(cs%isdw:cs%iedw,cs%jsdw-1:cs%jedw)) ; cs%OBCmask_v(:,:) = 0.0
6000 do j=g%jsd,g%jed ; do i=g%isd,g%ied
6001 cs%IareaT(i,j) = g%IareaT(i,j)
6002 cs%bathyT(i,j) = g%bathyT(i,j)
6003 cs%IareaT_OBCmask(i,j) = cs%IareaT(i,j)
6004 enddo ; enddo
6005
6006 ! Note: G%IdxCu & G%IdyCv may be valid for a smaller extent than CS%IdxCu & CS%IdyCv, even without
6007 ! wide halos.
6008 do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB
6009 cs%IdxCu(i,j) = g%IdxCu(i,j) ; cs%dy_Cu(i,j) = g%dy_Cu(i,j)
6010 cs%OBCmask_u(i,j) = g%OBCmaskCu(i,j)
6011 enddo ; enddo
6012 do j=g%JsdB,g%JedB ; do i=g%isd,g%ied
6013 cs%IdyCv(i,j) = g%IdyCv(i,j) ; cs%dx_Cv(i,j) = g%dx_Cv(i,j)
6014 cs%OBCmask_v(i,j) = g%OBCmaskCv(i,j)
6015 enddo ; enddo
6016
6017 ! This sets pressure force diagnostics on land, at coastlines and at OBC points to zero.
6018 if (mask_coastal_pressure_force) then
6019 do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB
6020 cs%IdxCu(i,j) = g%IdxCu_OBCmask(i,j)
6021 enddo ; enddo
6022 do j=g%JsdB,g%JedB ; do i=g%isd,g%ied
6023 cs%IdyCv(i,j) = g%IdyCv_OBCmask(i,j)
6024 enddo ; enddo
6025 endif
6026
6027 if (associated(obc)) then
6028 ! Set up information about the location and nature of the open boundary condition points.
6029 call initialize_bt_obc(obc, cs%BT_OBC, g, cs)
6030
6031 ! Update IareaT_OBCmask so that nothing changes outside of the OBC (problem for interior OBCs only)
6032 if (.not.cs%exterior_OBC_bug) then
6033 if (cs%BT_OBC%u_OBCs_on_PE) then
6034 do j=jsd,jed ; do i=isd,ied
6035 if (cs%BT_OBC%u_OBC_type(i-1,j) > 0) cs%IareaT_OBCmask(i,j) = 0.0 ! OBC_DIRECTION_E
6036 if (cs%BT_OBC%u_OBC_type(i,j) < 0) cs%IareaT_OBCmask(i,j) = 0.0 ! OBC_DIRECTION_W
6037 enddo ; enddo
6038 endif
6039 if (cs%BT_OBC%v_OBCs_on_PE) then
6040 do j=jsd,jed ; do i=isd,ied
6041 if (cs%BT_OBC%v_OBC_type(i,j-1) > 0) cs%IareaT_OBCmask(i,j) = 0.0 ! OBC_DIRECTION_N
6042 if (cs%BT_OBC%v_OBC_type(i,j) < 0) cs%IareaT_OBCmask(i,j) = 0.0 ! OBC_DIRECTION_S
6043 enddo ; enddo
6044 endif
6045 endif
6046
6047 ! Set masks to avoid changing velocities at OBC points.
6048 if (cs%BT_OBC%u_OBCs_on_PE) then
6049 do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB ; if (cs%BT_OBC%u_OBC_type(i,j) /= 0) then
6050 cs%OBCmask_u(i,j) = 0.0 ; cs%IdxCu(i,j) = 0.0
6051 endif ; enddo ; enddo
6052 endif
6053 if (cs%BT_OBC%v_OBCs_on_PE) then
6054 do j=g%JsdB,g%JedB ; do i=g%isd,g%ied ; if (cs%BT_OBC%v_OBC_type(i,j) /= 0) then
6055 cs%OBCmask_v(i,j) = 0.0 ; cs%IdyCv(i,j) = 0.0
6056 endif ; enddo ; enddo
6057 endif
6058
6059 cs%integral_OBCs = cs%integral_BT_cont .and. open_boundary_query(obc, apply_open_obc=.true.)
6060 else ! There are no OBC points anywhere.
6061 cs%BT_OBC%u_OBCs_on_PE = .false.
6062 cs%BT_OBC%v_OBCs_on_PE = .false.
6063 cs%integral_OBCs = .false.
6064 endif
6065
6066 call create_group_pass(pass_static_data, cs%IareaT, cs%BT_domain, to_all)
6067 call create_group_pass(pass_static_data, cs%bathyT, cs%BT_domain, to_all)
6068 call create_group_pass(pass_static_data, cs%IareaT_OBCmask, cs%BT_domain, to_all)
6069 call create_group_pass(pass_static_data, cs%IdxCu, cs%IdyCv, cs%BT_domain, to_all+scalar_pair)
6070 call create_group_pass(pass_static_data, cs%dy_Cu, cs%dx_Cv, cs%BT_domain, to_all+scalar_pair)
6071 call create_group_pass(pass_static_data, cs%OBCmask_u, cs%OBCmask_v, cs%BT_domain, to_all+scalar_pair)
6072 call do_group_pass(pass_static_data, cs%BT_domain)
6073
6074 ! Determine the weights to use for the thicknesses when calculating PV for use in the Coriolis terms
6075 allocate(cs%q_wt(4,cs%isdw-1:cs%iedw,cs%jsdw-1:cs%jedw), source=0.0)
6076 do j=js-1,je ; do i=is-1,ie
6077 if (g%mask2dT(i,j) + g%mask2dT(i,j+1) + g%mask2dT(i+1,j) + g%mask2dT(i+1,j+1) > 0.) then
6078 cs%q_wt(1,i,j) = g%areaT(i,j) ; cs%q_wt(2,i,j) = g%areaT(i+1,j)
6079 cs%q_wt(3,i,j) = g%areaT(i,j+1) ; cs%q_wt(4,i,j) = g%areaT(i+1,j+1)
6080 else
6081 cs%q_wt(1:4,i,j) = 0.0
6082 endif
6083 enddo ; enddo
6084
6085 if (cs%interior_OBC_PV .and. (cs%BT_OBC%u_OBCs_on_PE .or. cs%BT_OBC%v_OBCs_on_PE)) then
6086 ! Reset the potential vorticity at OBC vertices as a masked weighted average.
6087 do j=js-1,je ; do i=is-1,ie
6088 if ((g%mask2dT(i,j) + g%mask2dT(i,j+1) + g%mask2dT(i+1,j) + g%mask2dT(i+1,j+1) > 0.) .and. &
6089 ((abs(cs%BT_OBC%u_OBC_type(i,j)) > 0) .or. (abs(cs%BT_OBC%u_OBC_type(i,j+1)) > 0) .or. &
6090 (abs(cs%BT_OBC%v_OBC_type(i,j)) > 0) .or. (abs(cs%BT_OBC%v_OBC_type(i+1,j)) > 0)) ) then
6091 ! This is an OBC vertex, so use an area weighted masked average and avoid external values.
6092 cs%q_wt(1,i,j) = g%mask2dT(i,j) * g%areaT(i,j)
6093 cs%q_wt(2,i,j) = g%mask2dT(i+1,j) * g%areaT(i+1,j)
6094 cs%q_wt(3,i,j) = g%mask2dT(i,j+1) * g%areaT(i,j+1)
6095 cs%q_wt(4,i,j) = g%mask2dT(i+1,j+1) * g%areaT(i+1,j+1)
6096
6097 ! The following block is the equivalent of shifting weights inward across OBC points. With
6098 ! two OBCs in a line, it gives weights of about 1/2 and 1/2 to the interior points. At a
6099 ! peninsula-like corner between two OBCs it gives weights of about 3/8, 1/4 and 3/8 for the
6100 ! 3 interior points. At a bay-liek corner there is only one interior point with a weight of 1.
6101 ! The masking above zeros out the weights for exterior points.
6102 if (cs%BT_OBC%u_OBC_type(i,j) > 0) then ! Eastern OBC in the u-point to the south
6103 cs%q_wt(1,i,j) = cs%q_wt(1,i,j) + 0.5*g%mask2dT(i,j)*g%areaT(i,j) ! already CS%q_wt(2,I,J) = 0.0
6104 elseif (cs%BT_OBC%u_OBC_type(i,j) < 0) then ! Western OBC in the u-point to the south
6105 cs%q_wt(2,i,j) = cs%q_wt(2,i,j) + 0.5*g%mask2dT(i+1,j)*g%areaT(i+1,j) ! already CS%q_wt(1,I,J) = 0.0
6106 endif
6107 if (cs%BT_OBC%u_OBC_type(i,j+1) > 0) then ! Eastern OBC in the u-point to the north
6108 cs%q_wt(3,i,j) = cs%q_wt(3,i,j) + 0.5*g%mask2dT(i,j+1)*g%areaT(i,j+1) ! already CS%q_wt(4,I,J) = 0.0
6109 elseif (cs%BT_OBC%u_OBC_type(i,j+1) < 0) then ! Western OBC in the u-point to the north
6110 cs%q_wt(4,i,j) = cs%q_wt(4,i,j) + 0.5*g%mask2dT(i+1,j+1)*g%areaT(i+1,j+1) ! already CS%q_wt(3,I,J) = 0.0
6111 endif
6112 if (cs%BT_OBC%v_OBC_type(i,j) > 0) then ! Northern OBC in the v-point to the west
6113 cs%q_wt(1,i,j) = cs%q_wt(1,i,j) + 0.5*g%mask2dT(i,j)*g%areaT(i,j) ! already CS%q_wt(3,I,J) = 0.0
6114 elseif (cs%BT_OBC%v_OBC_type(i,j) < 0) then ! Southern OBC in the v-point to the west
6115 cs%q_wt(3,i,j) = cs%q_wt(3,i,j) + 0.5*g%mask2dT(i,j+1)*g%areaT(i,j+1) ! already CS%q_wt(1,I,J) = 0.0
6116 endif
6117 if (cs%BT_OBC%v_OBC_type(i+1,j) > 0) then ! Northern OBC in the v-point to the west
6118 cs%q_wt(2,i,j) = cs%q_wt(2,i,j) + 0.5*g%mask2dT(i+1,j)*g%areaT(i+1,j) ! already CS%q_wt(4,I,J) = 0.0
6119 elseif (cs%BT_OBC%v_OBC_type(i+1,j) < 0) then ! Southern OBC in the v-point to the west
6120 cs%q_wt(4,i,j) = cs%q_wt(4,i,j) + 0.5*g%mask2dT(i+1,j+1)*g%areaT(i+1,j+1) ! already CS%q_wt(2,I,J) = 0.0
6121 endif
6122 endif
6123 enddo ; enddo
6124 endif
6125
6126 if (cs%linearized_BT_PV) then
6127 allocate(cs%q_D(cs%isdw-1:cs%iedw,cs%jsdw-1:cs%jedw), source=0.0)
6128 allocate(cs%D_u_Cor(cs%isdw-1:cs%iedw,cs%jsdw:cs%jedw), source=0.0)
6129 allocate(cs%D_v_Cor(cs%isdw:cs%iedw,cs%jsdw-1:cs%jedw), source=0.0)
6130
6131 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
6132
6133 do j=js,je ; do i=is-1,ie
6134 cs%D_u_Cor(i,j) = 0.5 * ( max(g%meanSL(i+1,j) + g%bathyT(i+1,j), 0.0) &
6135 + max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) ) * z_to_h
6136 enddo ; enddo
6137 if (cs%interior_OBC_PV .and. cs%BT_OBC%u_OBCs_on_PE) then ; do j=js,je ; do i=is-1,ie
6138 if (cs%BT_OBC%u_OBC_type(i,j) < 0) & ! Western boundary condition
6139 cs%D_u_Cor(i,j) = max(g%meanSL(i+1,j) + g%bathyT(i+1,j), 0.0) * z_to_h
6140 if (cs%BT_OBC%u_OBC_type(i,j) > 0) & ! Eastern boundary condition
6141 cs%D_u_Cor(i,j) = max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) * z_to_h
6142 enddo ; enddo ; endif
6143
6144 do j=js-1,je ; do i=is,ie
6145 cs%D_v_Cor(i,j) = 0.5 * ( max(g%meanSL(i,j+1) + g%bathyT(i,j+1), 0.0) &
6146 + max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) ) * z_to_h
6147 enddo ; enddo
6148 if (cs%interior_OBC_PV .and. cs%BT_OBC%v_OBCs_on_PE) then ; do j=js-1,je ; do i=is,ie
6149 if (cs%BT_OBC%v_OBC_type(i,j) < 0) & ! Southern boundary condition
6150 cs%D_v_Cor(i,j) = max(g%meanSL(i,j+1) + g%bathyT(i,j+1), 0.0) * z_to_h
6151 if (cs%BT_OBC%v_OBC_type(i,j) > 0) & ! Northern boundary condition
6152 cs%D_v_Cor(i,j) = max(g%meanSL(i,j) + g%bathyT(i,j), 0.0) * z_to_h
6153 enddo ; enddo ; endif
6154
6155 h_a_neglect = gv%H_subroundoff * 1.0 * us%m_to_L**2
6156 do j=js-1,je ; do i=is-1,ie
6157 if ((cs%q_wt(1,i,j) + cs%q_wt(4,i,j)) + (cs%q_wt(2,i,j) + cs%q_wt(3,i,j)) > 0.) then
6158 cs%q_D(i,j) = 0.25 * (cs%BT_Coriolis_scale * g%CoriolisBu(i,j)) * &
6159 ((cs%q_wt(1,i,j) + cs%q_wt(4,i,j)) + (cs%q_wt(2,i,j) + cs%q_wt(3,i,j))) / &
6160 max(z_to_h * (((cs%q_wt(1,i,j) * max(g%meanSL(i,j) + g%bathyT(i,j), 0.0)) + &
6161 (cs%q_wt(4,i,j) * max(g%meanSL(i+1,j+1) + g%bathyT(i+1,j+1), 0.0))) + &
6162 ((cs%q_wt(2,i,j) * max(g%meanSL(i+1,j) + g%bathyT(i+1,j), 0.0)) + &
6163 (cs%q_wt(3,i,j) * max(g%meanSL(i,j+1) + g%bathyT(i,j+1), 0.0)))), &
6164 h_a_neglect)
6165 else ! All four h points are masked out so q_D(I,J) is meaningless
6166 cs%q_D(i,j) = 0.
6167 endif
6168 enddo ; enddo
6169
6170 ! With very wide halos, q and D need to be calculated on the available data
6171 ! domain and then updated onto the full computational domain.
6172 call create_group_pass(pass_q_d_cor, cs%q_D, cs%BT_Domain, to_all, position=corner)
6173 call create_group_pass(pass_q_d_cor, cs%D_u_Cor, cs%D_v_Cor, cs%BT_Domain, &
6174 to_all+scalar_pair)
6175 call do_group_pass(pass_q_d_cor, cs%BT_Domain)
6176 endif
6177
6178 if (cs%linear_wave_drag) then
6179 allocate(cs%lin_drag_u(isdb:iedb,jsd:jed), source=0.0)
6180 allocate(cs%lin_drag_v(isd:ied,jsdb:jedb), source=0.0)
6181
6182 if (len_trim(wave_drag_file) > 0) then
6183 inputdir = "." ; call get_param(param_file, mdl, "INPUTDIR", inputdir)
6184 wave_drag_file = trim(slasher(inputdir))//trim(wave_drag_file)
6185 call log_param(param_file, mdl, "INPUTDIR/BT_WAVE_DRAG_FILE", wave_drag_file)
6186
6187 if (len_trim(wave_drag_u) > 0 .and. len_trim(wave_drag_v) > 0) then
6188 call mom_read_data(wave_drag_file, wave_drag_u, cs%lin_drag_u, g%Domain, &
6189 position=east_face, scale=wave_drag_scale*gv%m_to_H*us%T_to_s)
6190 call mom_read_data(wave_drag_file, wave_drag_v, cs%lin_drag_v, g%Domain, &
6191 position=north_face, scale=wave_drag_scale*gv%m_to_H*us%T_to_s)
6192 call pass_vector(cs%lin_drag_u, cs%lin_drag_v, g%domain, direction=to_all+scalar_pair)
6193 else
6194 allocate(lin_drag_h(isd:ied,jsd:jed), source=0.0)
6195
6196 call mom_read_data(wave_drag_file, wave_drag_var, lin_drag_h, g%Domain, scale=gv%m_to_H*us%T_to_s)
6197 call pass_var(lin_drag_h, g%Domain)
6198 do j=js,je ; do i=is-1,ie
6199 cs%lin_drag_u(i,j) = wave_drag_scale * 0.5 * (lin_drag_h(i,j) + lin_drag_h(i+1,j))
6200 enddo ; enddo
6201 do j=js-1,je ; do i=is,ie
6202 cs%lin_drag_v(i,j) = wave_drag_scale * 0.5 * (lin_drag_h(i,j) + lin_drag_h(i,j+1))
6203 enddo ; enddo
6204 deallocate(lin_drag_h)
6205 endif ! len_trim(wave_drag_u) > 0 .and. len_trim(wave_drag_v) > 0
6206 endif ! len_trim(wave_drag_file) > 0
6207 endif ! CS%linear_wave_drag
6208
6209 ! Initialize streaming band-pass filters and frequency-dependent drag
6210 if (cs%use_filter) then
6211 call filt_init(param_file, us, cs%Filt_CS_u, restart_cs)
6212 call filt_init(param_file, us, cs%Filt_CS_v, restart_cs)
6213 endif
6214
6215 if (cs%use_filter .and. cs%linear_freq_drag) then
6216 if (.not.cs%linear_wave_drag .and. len_trim(wave_drag_file) > 0) then
6217 inputdir = "." ; call get_param(param_file, mdl, "INPUTDIR", inputdir)
6218 wave_drag_file = trim(slasher(inputdir))//trim(wave_drag_file)
6219 call log_param(param_file, mdl, "INPUTDIR/BT_WAVE_DRAG_FILE", wave_drag_file)
6220 endif
6221 call wave_drag_init(param_file, wave_drag_file, g, gv, us, cs%Drag_CS)
6222 endif
6223
6224 cs%dtbt_fraction = 0.98 ; if (dtbt_input < 0.0) cs%dtbt_fraction = -dtbt_input
6225
6226 dtbt_restart = -1.0
6227 if (query_initialized(cs%dtbt, "DTBT", restart_cs)) then
6228 dtbt_restart = cs%dtbt
6229 endif
6230
6231 ! Estimate the maximum stable barotropic time step.
6232 gtot_estimate = 0.0
6233 if (gv%Boussinesq) then
6234 do k=1,gv%ke ; gtot_estimate = gtot_estimate + gv%H_to_Z*gv%g_prime(k) ; enddo
6235 else
6236 h_to_z = gv%H_to_RZ / cs%Rho_BT_lin
6237 do k=1,gv%ke ; gtot_estimate = gtot_estimate + h_to_z*gv%g_prime(k) ; enddo
6238 endif
6239
6240 ! CS%dtbt calculated here by set_dtbt is only used when dtbt is not reset during the run, i.e. DTBT_RESET_PERIOD<0.
6241 call set_dtbt(g, gv, us, cs, gtot_est=gtot_estimate, ssh_add=ssh_extra)
6242
6243 if (dtbt_input > 0.0) then
6244 cs%dtbt = us%s_to_T * dtbt_input
6245 elseif (dtbt_restart > 0.0) then
6246 cs%dtbt = dtbt_restart
6247 endif
6248
6249 calc_dtbt = .true. ; if ((dtbt_restart > 0.0) .and. (dtbt_input > 0.0)) calc_dtbt = .false.
6250
6251 call log_param(param_file, mdl, "DTBT as used", cs%dtbt, units="s", unscale=us%T_to_s)
6252 call log_param(param_file, mdl, "estimated maximum DTBT", cs%dtbt_max, units="s", unscale=us%T_to_s)
6253
6254 ! ubtav and vbtav, and perhaps ubt_IC and vbt_IC, are allocated and
6255 ! initialized in register_barotropic_restarts.
6256
6257 if (gv%Boussinesq) then
6258 thickness_units = "m" ; flux_units = "m3 s-1"
6259 else
6260 thickness_units = "kg m-2" ; flux_units = "kg s-1"
6261 endif
6262
6263 cs%id_PFu_bt = register_diag_field('ocean_model', 'PFuBT', diag%axesCu1, time, &
6264 'Zonal Anomalous Barotropic Pressure Force Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6265 cs%id_PFv_bt = register_diag_field('ocean_model', 'PFvBT', diag%axesCv1, time, &
6266 'Meridional Anomalous Barotropic Pressure Force Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6267 cs%id_etaPF_anom = register_diag_field('ocean_model', 'etaPF_anom', diag%axesT1, time, &
6268 'Eta anomalies used for pressure forces averaged over a baroclinic timestep', &
6269 thickness_units, conversion=gv%H_to_MKS)
6270 if (cs%linear_wave_drag .or. (cs%use_filter .and. cs%linear_freq_drag)) then
6271 cs%id_LDu_bt = register_diag_field('ocean_model', 'WaveDraguBT', diag%axesCu1, time, &
6272 'Zonal Barotropic Linear Wave Drag Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6273 cs%id_LDv_bt = register_diag_field('ocean_model', 'WaveDragvBT', diag%axesCv1, time, &
6274 'Meridional Barotropic Linear Wave Drag Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6275 endif
6276 cs%id_Coru_bt = register_diag_field('ocean_model', 'CoruBT', diag%axesCu1, time, &
6277 'Zonal Barotropic Coriolis Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6278 cs%id_Corv_bt = register_diag_field('ocean_model', 'CorvBT', diag%axesCv1, time, &
6279 'Meridional Barotropic Coriolis Acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6280 cs%id_uaccel = register_diag_field('ocean_model', 'u_accel_bt', diag%axesCu1, time, &
6281 'Barotropic zonal acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6282 cs%id_vaccel = register_diag_field('ocean_model', 'v_accel_bt', diag%axesCv1, time, &
6283 'Barotropic meridional acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6284 cs%id_ubtforce = register_diag_field('ocean_model', 'ubtforce', diag%axesCu1, time, &
6285 'Barotropic zonal acceleration from baroclinic terms', 'm s-2', conversion=us%L_T2_to_m_s2)
6286 cs%id_vbtforce = register_diag_field('ocean_model', 'vbtforce', diag%axesCv1, time, &
6287 'Barotropic meridional acceleration from baroclinic terms', 'm s-2', conversion=us%L_T2_to_m_s2)
6288 cs%id_ubtdt = register_diag_field('ocean_model', 'ubt_dt', diag%axesCu1, time, &
6289 'Barotropic zonal acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6290 cs%id_vbtdt = register_diag_field('ocean_model', 'vbt_dt', diag%axesCv1, time, &
6291 'Barotropic meridional acceleration', 'm s-2', conversion=us%L_T2_to_m_s2)
6292
6293 cs%id_eta_bt = register_diag_field('ocean_model', 'eta_bt', diag%axesT1, time, &
6294 'Barotropic end SSH', thickness_units, conversion=gv%H_to_MKS)
6295 cs%id_ubt = register_diag_field('ocean_model', 'ubt', diag%axesCu1, time, &
6296 'Barotropic end zonal velocity', 'm s-1', conversion=us%L_T_to_m_s)
6297 cs%id_vbt = register_diag_field('ocean_model', 'vbt', diag%axesCv1, time, &
6298 'Barotropic end meridional velocity', 'm s-1', conversion=us%L_T_to_m_s)
6299 cs%id_eta_st = register_diag_field('ocean_model', 'eta_st', diag%axesT1, time, &
6300 'Barotropic start SSH', thickness_units, conversion=gv%H_to_MKS)
6301 cs%id_ubt_st = register_diag_field('ocean_model', 'ubt_st', diag%axesCu1, time, &
6302 'Barotropic start zonal velocity', 'm s-1', conversion=us%L_T_to_m_s)
6303 cs%id_vbt_st = register_diag_field('ocean_model', 'vbt_st', diag%axesCv1, time, &
6304 'Barotropic start meridional velocity', 'm s-1', conversion=us%L_T_to_m_s)
6305 cs%id_ubtav = register_diag_field('ocean_model', 'ubtav', diag%axesCu1, time, &
6306 'Barotropic time-average zonal velocity', 'm s-1', conversion=us%L_T_to_m_s)
6307 cs%id_vbtav = register_diag_field('ocean_model', 'vbtav', diag%axesCv1, time, &
6308 'Barotropic time-average meridional velocity', 'm s-1', conversion=us%L_T_to_m_s)
6309 cs%id_eta_cor = register_diag_field('ocean_model', 'eta_cor', diag%axesT1, time, &
6310 'Corrective mass or volume flux within a timestep', thickness_units, conversion=gv%H_to_MKS)
6311 cs%id_visc_rem_u = register_diag_field('ocean_model', 'visc_rem_u', diag%axesCuL, time, &
6312 'Viscous remnant at u', 'nondim')
6313 cs%id_visc_rem_v = register_diag_field('ocean_model', 'visc_rem_v', diag%axesCvL, time, &
6314 'Viscous remnant at v', 'nondim')
6315 cs%id_bt_rem_u = register_diag_field('ocean_model', 'bt_rem_u', diag%axesCu1, time, &
6316 'Barotropic viscous remnant per barotropic step at u', 'nondim')
6317 cs%id_bt_rem_v = register_diag_field('ocean_model', 'bt_rem_v', diag%axesCv1, time, &
6318 'Barotropic viscous remnant per barotropic step at v', 'nondim')
6319 cs%id_gtotn = register_diag_field('ocean_model', 'gtot_n', diag%axesT1, time, &
6320 'gtot to North', 'm s-2', conversion=gv%m_to_H*(us%L_T_to_m_s**2))
6321 cs%id_gtots = register_diag_field('ocean_model', 'gtot_s', diag%axesT1, time, &
6322 'gtot to South', 'm s-2', conversion=gv%m_to_H*(us%L_T_to_m_s**2))
6323 cs%id_gtote = register_diag_field('ocean_model', 'gtot_e', diag%axesT1, time, &
6324 'gtot to East', 'm s-2', conversion=gv%m_to_H*(us%L_T_to_m_s**2))
6325 cs%id_gtotw = register_diag_field('ocean_model', 'gtot_w', diag%axesT1, time, &
6326 'gtot to West', 'm s-2', conversion=gv%m_to_H*(us%L_T_to_m_s**2))
6327 cs%id_eta_hifreq = register_diag_field('ocean_model', 'eta_hifreq', diag%axesT1, time, &
6328 'High Frequency Barotropic SSH', thickness_units, conversion=gv%H_to_MKS)
6329 cs%id_ubt_hifreq = register_diag_field('ocean_model', 'ubt_hifreq', diag%axesCu1, time, &
6330 'High Frequency Barotropic zonal velocity', 'm s-1', conversion=us%L_T_to_m_s)
6331 cs%id_vbt_hifreq = register_diag_field('ocean_model', 'vbt_hifreq', diag%axesCv1, time, &
6332 'High Frequency Barotropic meridional velocity', 'm s-1', conversion=us%L_T_to_m_s)
6333 ! if (.not.CS%BT_project_velocity) & ! The following diagnostic is redundant with BT_project_velocity.
6334 cs%id_eta_pred_hifreq = register_diag_field('ocean_model', 'eta_pred_hifreq', diag%axesT1, time, &
6335 'High Frequency Predictor Barotropic SSH', thickness_units, conversion=gv%H_to_MKS)
6336 cs%id_etaPF_hifreq = register_diag_field('ocean_model', 'etaPF_hifreq', diag%axesT1, time, &
6337 'High Frequency Barotropic SSH anomalies used for pressure forces', thickness_units, conversion=gv%H_to_MKS)
6338 cs%id_uhbt_hifreq = register_diag_field('ocean_model', 'uhbt_hifreq', diag%axesCu1, time, &
6339 'High Frequency Barotropic zonal transport', &
6340 'm3 s-1', conversion=gv%H_to_m*us%L_to_m*us%L_T_to_m_s)
6341 cs%id_vhbt_hifreq = register_diag_field('ocean_model', 'vhbt_hifreq', diag%axesCv1, time, &
6342 'High Frequency Barotropic meridional transport', &
6343 'm3 s-1', conversion=gv%H_to_m*us%L_to_m*us%L_T_to_m_s)
6344 cs%id_frhatu = register_diag_field('ocean_model', 'frhatu', diag%axesCuL, time, &
6345 'Fractional thickness of layers in u-columns', 'nondim')
6346 cs%id_frhatv = register_diag_field('ocean_model', 'frhatv', diag%axesCvL, time, &
6347 'Fractional thickness of layers in v-columns', 'nondim')
6348 cs%id_frhatu1 = register_diag_field('ocean_model', 'frhatu1', diag%axesCuL, time, &
6349 'Predictor Fractional thickness of layers in u-columns', 'nondim')
6350 cs%id_frhatv1 = register_diag_field('ocean_model', 'frhatv1', diag%axesCvL, time, &
6351 'Predictor Fractional thickness of layers in v-columns', 'nondim')
6352 cs%id_uhbt = register_diag_field('ocean_model', 'uhbt', diag%axesCu1, time, &
6353 'Barotropic zonal transport averaged over a baroclinic step', &
6354 'm3 s-1', conversion=gv%H_to_m*us%L_to_m*us%L_T_to_m_s)
6355 cs%id_vhbt = register_diag_field('ocean_model', 'vhbt', diag%axesCv1, time, &
6356 'Barotropic meridional transport averaged over a baroclinic step', &
6357 'm3 s-1', conversion=gv%H_to_m*us%L_to_m*us%L_T_to_m_s)
6358
6359 if (use_bt_cont_type) then
6360 cs%id_BTC_FA_u_EE = register_diag_field('ocean_model', 'BTC_FA_u_EE', diag%axesCu1, time, &
6361 'BTCont type far east face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6362 cs%id_BTC_FA_u_E0 = register_diag_field('ocean_model', 'BTC_FA_u_E0', diag%axesCu1, time, &
6363 'BTCont type near east face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6364 cs%id_BTC_FA_u_WW = register_diag_field('ocean_model', 'BTC_FA_u_WW', diag%axesCu1, time, &
6365 'BTCont type far west face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6366 cs%id_BTC_FA_u_W0 = register_diag_field('ocean_model', 'BTC_FA_u_W0', diag%axesCu1, time, &
6367 'BTCont type near west face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6368 cs%id_BTC_ubt_EE = register_diag_field('ocean_model', 'BTC_ubt_EE', diag%axesCu1, time, &
6369 'BTCont type far east velocity', 'm s-1', conversion=us%L_T_to_m_s)
6370 cs%id_BTC_ubt_WW = register_diag_field('ocean_model', 'BTC_ubt_WW', diag%axesCu1, time, &
6371 'BTCont type far west velocity', 'm s-1', conversion=us%L_T_to_m_s)
6372 ! This is a specialized diagnostic that is not being made widely available (yet).
6373 ! CS%id_BTC_FA_u_rat0 = register_diag_field('ocean_model', 'BTC_FA_u_rat0', diag%axesCu1, Time, &
6374 ! 'BTCont type ratio of near east and west face areas', 'nondim')
6375 cs%id_BTC_FA_v_NN = register_diag_field('ocean_model', 'BTC_FA_v_NN', diag%axesCv1, time, &
6376 'BTCont type far north face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6377 cs%id_BTC_FA_v_N0 = register_diag_field('ocean_model', 'BTC_FA_v_N0', diag%axesCv1, time, &
6378 'BTCont type near north face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6379 cs%id_BTC_FA_v_SS = register_diag_field('ocean_model', 'BTC_FA_v_SS', diag%axesCv1, time, &
6380 'BTCont type far south face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6381 cs%id_BTC_FA_v_S0 = register_diag_field('ocean_model', 'BTC_FA_v_S0', diag%axesCv1, time, &
6382 'BTCont type near south face area', 'm2', conversion=us%L_to_m*gv%H_to_m)
6383 cs%id_BTC_vbt_NN = register_diag_field('ocean_model', 'BTC_vbt_NN', diag%axesCv1, time, &
6384 'BTCont type far north velocity', 'm s-1', conversion=us%L_T_to_m_s)
6385 cs%id_BTC_vbt_SS = register_diag_field('ocean_model', 'BTC_vbt_SS', diag%axesCv1, time, &
6386 'BTCont type far south velocity', 'm s-1', conversion=us%L_T_to_m_s)
6387 ! This is a specialized diagnostic that is not being made widely available (yet).
6388 ! CS%id_BTC_FA_v_rat0 = register_diag_field('ocean_model', 'BTC_FA_v_rat0', diag%axesCv1, Time, &
6389 ! 'BTCont type ratio of near north and south face areas', 'nondim')
6390 ! CS%id_BTC_FA_h_rat0 = register_diag_field('ocean_model', 'BTC_FA_h_rat0', diag%axesT1, Time, &
6391 ! 'BTCont type maximum ratios of near face areas around cells', 'nondim')
6392 endif
6393 cs%id_uhbt0 = register_diag_field('ocean_model', 'uhbt0', diag%axesCu1, time, &
6394 'Barotropic zonal transport difference', 'm3 s-1', conversion=gv%H_to_m*us%L_to_m**2*us%s_to_T)
6395 cs%id_vhbt0 = register_diag_field('ocean_model', 'vhbt0', diag%axesCv1, time, &
6396 'Barotropic meridional transport difference', 'm3 s-1', conversion=gv%H_to_m*us%L_to_m**2*us%s_to_T)
6397 if (associated(obc)) then
6398 if (obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally) then
6399 cs%id_SSH_u_OBC = register_diag_field('ocean_model', 'SSH_u_OBC', diag%axesCu1, time, &
6400 'Outer sea surface height at u OBC points', 'm', conversion=us%Z_to_m)
6401 cs%id_SSH_v_OBC = register_diag_field('ocean_model', 'SSH_v_OBC', diag%axesCv1, time, &
6402 'Outer sea surface height at v OBC points', 'm', conversion=us%Z_to_m)
6403 cs%id_ubt_OBC = register_diag_field('ocean_model', 'ubt_OBC', diag%axesCu1, time, &
6404 'Outer u velocity at OBC points', 'm', conversion=us%L_T_to_m_s)
6405 cs%id_vbt_OBC = register_diag_field('ocean_model', 'vbt_OBC', diag%axesCv1, time, &
6406 'Outer v velocity at OBC points', 'm', conversion=us%L_T_to_m_s)
6407 endif
6408 endif
6409
6410 if (cs%id_frhatu1 > 0) allocate(cs%frhatu1(isdb:iedb,jsd:jed,nz), source=0.)
6411 if (cs%id_frhatv1 > 0) allocate(cs%frhatv1(isd:ied,jsdb:jedb,nz), source=0.)
6412
6413 if (.NOT.query_initialized(cs%ubtav,"ubtav",restart_cs) .or. &
6414 .NOT.query_initialized(cs%vbtav,"vbtav",restart_cs)) then
6415 call btcalc(h, g, gv, cs, may_use_default=.true.)
6416 cs%ubtav(:,:) = 0.0 ; cs%vbtav(:,:) = 0.0
6417 do k=1,nz ; do j=js,je ; do i=is-1,ie
6418 cs%ubtav(i,j) = cs%ubtav(i,j) + cs%frhatu(i,j,k) * u(i,j,k)
6419 enddo ; enddo ; enddo
6420 do k=1,nz ; do j=js-1,je ; do i=is,ie
6421 cs%vbtav(i,j) = cs%vbtav(i,j) + cs%frhatv(i,j,k) * v(i,j,k)
6422 enddo ; enddo ; enddo
6423 endif
6424
6425 if (cs%gradual_BT_ICs) then
6426 if (.NOT.query_initialized(cs%ubt_IC,"ubt_IC",restart_cs) .or. &
6427 .NOT.query_initialized(cs%vbt_IC,"vbt_IC",restart_cs)) then
6428 do j=js,je ; do i=is-1,ie ; cs%ubt_IC(i,j) = cs%ubtav(i,j) ; enddo ; enddo
6429 do j=js-1,je ; do i=is,ie ; cs%vbt_IC(i,j) = cs%vbtav(i,j) ; enddo ; enddo
6430 endif
6431 endif
6432 ! Calculate other constants which are used for btstep.
6433
6434 if (.not.cs%nonlin_stress) then
6435 z_to_h = gv%Z_to_H ; if (.not.gv%Boussinesq) z_to_h = gv%RZ_to_H * cs%Rho_BT_lin
6436 do j=js,je ; do i=is-1,ie
6437 htot = max(g%meanSL(i+1,j) + g%bathyT(i+1,j), 0.0) + max(g%meanSL(i,j) + g%bathyT(i,j), 0.0)
6438 if (g%OBCmaskCu(i,j) * htot > 0.) then
6439 cs%IDatu(i,j) = g%OBCmaskCu(i,j) * 2.0 / (z_to_h * htot)
6440 else ! Both neighboring H points are masked out or this is an OBC face so IDatu(I,j) is unused
6441 cs%IDatu(i,j) = 0.
6442 endif
6443 enddo ; enddo
6444 do j=js-1,je ; do i=is,ie
6445 htot = max(g%meanSL(i,j+1) + g%bathyT(i,j+1), 0.0) + max(g%meanSL(i,j) + g%bathyT(i,j), 0.0)
6446 if (g%OBCmaskCv(i,j) * htot > 0.) then
6447 cs%IDatv(i,j) = g%OBCmaskCv(i,j) * 2.0 / (z_to_h * htot)
6448 else ! Both neighboring H points are masked out or this is an OBC face so IDatv(i,J) is unused
6449 cs%IDatv(i,j) = 0.
6450 endif
6451 enddo ; enddo
6452 endif
6453
6454 call find_face_areas(datu, datv, g, gv, us, cs, ms, 1)
6455 if ((cs%bound_BT_corr) .and. .not.(use_bt_cont_type .and. cs%BT_cont_bounds)) then
6456 ! This is not used in most test cases. Were it ever to become more widely used, consider
6457 ! replacing maxvel with min(G%dxT(i,j),G%dyT(i,j)) * (CS%maxCFL_BT_cont*Idt) .
6458 do j=js,je ; do i=is,ie
6459 cs%eta_cor_bound(i,j) = g%IareaT(i,j) * 0.1 * cs%maxvel * &
6460 ((datu(i-1,j) + datu(i,j)) + (datv(i,j) + datv(i,j-1)))
6461 enddo ; enddo
6462 endif
6463
6464 if (cs%gradual_BT_ICs) &
6465 call create_group_pass(pass_bt_hbt_btav, cs%ubt_IC, cs%vbt_IC, g%Domain)
6466 call create_group_pass(pass_bt_hbt_btav, cs%ubtav, cs%vbtav, g%Domain)
6467 call do_group_pass(pass_bt_hbt_btav, g%Domain)
6468
6469 ! id_clock_pass = cpu_clock_id('(Ocean BT halo updates)', grain=CLOCK_ROUTINE)
6470 id_clock_calc_pre = cpu_clock_id('(Ocean BT pre-calcs only)', grain=clock_routine)
6471 id_clock_pass_pre = cpu_clock_id('(Ocean BT pre-step halo updates)', grain=clock_routine)
6472 id_clock_calc = cpu_clock_id('(Ocean BT stepping calcs only)', grain=clock_routine)
6473 id_clock_pass_step = cpu_clock_id('(Ocean BT stepping halo updates)', grain=clock_routine)
6474 id_clock_calc_post = cpu_clock_id('(Ocean BT post-calcs only)', grain=clock_routine)
6475 id_clock_pass_post = cpu_clock_id('(Ocean BT post-step halo updates)', grain=clock_routine)
6476 if (dtbt_input <= 0.0) &
6477 id_clock_sync = cpu_clock_id('(Ocean BT global synch)', grain=clock_routine)
6478
6479end subroutine barotropic_init
6480
6481!> Copies ubtav and vbtav from private type into arrays
6482subroutine barotropic_get_tav(CS, ubtav, vbtav, G, US)
6483 type(barotropic_cs), intent(in) :: cs !< Barotropic control structure
6484 type(ocean_grid_type), intent(in) :: g !< Grid structure
6485 real, dimension(SZIB_(G),SZJ_(G)), intent(inout) :: ubtav !< Zonal barotropic velocity averaged
6486 !! over a baroclinic timestep [L T-1 ~> m s-1]
6487 real, dimension(SZI_(G),SZJB_(G)), intent(inout) :: vbtav !< Meridional barotropic velocity averaged
6488 !! over a baroclinic timestep [L T-1 ~> m s-1]
6489 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
6490 ! Local variables
6491 integer :: i,j
6492
6493 do j=g%jsc,g%jec ; do i=g%isc-1,g%iec
6494 ubtav(i,j) = cs%ubtav(i,j)
6495 enddo ; enddo
6496
6497 do j=g%jsc-1,g%jec ; do i=g%isc,g%iec
6498 vbtav(i,j) = cs%vbtav(i,j)
6499 enddo ; enddo
6500
6501end subroutine barotropic_get_tav
6502
6503
6504!> Clean up the barotropic control structure.
6505subroutine barotropic_end(CS)
6506 type(barotropic_cs), intent(inout) :: cs !< Control structure to clear out.
6507
6508 call destroy_bt_obc(cs%BT_OBC)
6509
6510 ! Allocated in barotropic_init, called in timestep initialization
6511 dealloc_(cs%ua_polarity) ; dealloc_(cs%va_polarity)
6512 dealloc_(cs%IDatu) ; dealloc_(cs%IDatv)
6513 if (allocated(cs%eta_cor_bound)) deallocate(cs%eta_cor_bound)
6514 dealloc_(cs%eta_cor)
6515 dealloc_(cs%bathyT) ; dealloc_(cs%IareaT)
6516 dealloc_(cs%frhatu) ; dealloc_(cs%frhatv)
6517 dealloc_(cs%OBCmask_u) ; dealloc_(cs%OBCmask_v)
6518 dealloc_(cs%IdxCu) ; dealloc_(cs%IdyCv)
6519 dealloc_(cs%dy_Cu) ; dealloc_(cs%dx_Cv)
6520
6521 if (allocated(cs%frhatu1)) deallocate(cs%frhatu1)
6522 if (allocated(cs%frhatv1)) deallocate(cs%frhatv1)
6523 if (allocated(cs%IareaT_OBCmask)) deallocate(cs%IareaT_OBCmask)
6524
6525 if (allocated(cs%q_D)) deallocate(cs%q_D)
6526 if (allocated(cs%D_u_Cor)) deallocate(cs%D_u_Cor)
6527 if (allocated(cs%D_v_Cor)) deallocate(cs%D_v_Cor)
6528 if (allocated(cs%ubt_IC)) deallocate(cs%ubt_IC)
6529 if (allocated(cs%vbt_IC)) deallocate(cs%vbt_IC)
6530 if (allocated(cs%lin_drag_u)) deallocate(cs%lin_drag_u)
6531 if (allocated(cs%lin_drag_v)) deallocate(cs%lin_drag_v)
6532
6533 if (associated(cs%debug_BT_HI)) deallocate(cs%debug_BT_HI)
6534 call deallocate_mom_domain(cs%BT_domain)
6535
6536 ! Allocated in restart registration, prior to timestep initialization
6537 dealloc_(cs%ubtav) ; dealloc_(cs%vbtav)
6538end subroutine barotropic_end
6539
6540!> This subroutine is used to register any fields from MOM_barotropic.F90
6541!! that should be written to or read from the restart file.
6542subroutine register_barotropic_restarts(HI, GV, US, param_file, CS, restart_CS)
6543 type(hor_index_type), intent(in) :: hi !< A horizontal index type structure.
6544 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
6545 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
6546 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
6547 type(barotropic_cs), intent(inout) :: cs !< Barotropic control structure
6548 type(mom_restart_cs), intent(inout) :: restart_cs !< MOM restart control structure
6549
6550 ! Local variables
6551 type(vardesc) :: vd(3)
6552 character(len=40) :: mdl = "MOM_barotropic" ! This module's name.
6553 integer :: n_filters !< Number of streaming band-pass filters to be used in the simulation.
6554 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
6555
6556 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed
6557 isdb = hi%IsdB ; iedb = hi%IedB ; jsdb = hi%JsdB ; jedb = hi%JedB
6558
6559 call get_param(param_file, mdl, "GRADUAL_BT_ICS", cs%gradual_BT_ICs, &
6560 "If true, adjust the initial conditions for the "//&
6561 "barotropic solver to the values from the layered "//&
6562 "solution over a whole timestep instead of instantly. "//&
6563 "This is a decent approximation to the inclusion of "//&
6564 "sum(u dh_dt) while also correcting for truncation errors.", &
6565 default=.false., do_not_log=.true.)
6566
6567 alloc_(cs%ubtav(isdb:iedb,jsd:jed)) ; cs%ubtav(:,:) = 0.0
6568 alloc_(cs%vbtav(isd:ied,jsdb:jedb)) ; cs%vbtav(:,:) = 0.0
6569 if (cs%gradual_BT_ICs) then
6570 allocate(cs%ubt_IC(isdb:iedb,jsd:jed), source=0.0)
6571 allocate(cs%vbt_IC(isd:ied,jsdb:jedb), source=0.0)
6572 endif
6573
6574 vd(2) = var_desc("ubtav","m s-1","Time mean barotropic zonal velocity", &
6575 hor_grid='u', z_grid='1')
6576 vd(3) = var_desc("vbtav","m s-1","Time mean barotropic meridional velocity",&
6577 hor_grid='v', z_grid='1')
6578 call register_restart_pair(cs%ubtav, cs%vbtav, vd(2), vd(3), .false., restart_cs, &
6579 conversion=us%L_T_to_m_s)
6580
6581 if (cs%gradual_BT_ICs) then
6582 vd(2) = var_desc("ubt_IC", "m s-1", &
6583 longname="Next initial condition for the barotropic zonal velocity", &
6584 hor_grid='u', z_grid='1')
6585 vd(3) = var_desc("vbt_IC", "m s-1", &
6586 longname="Next initial condition for the barotropic meridional velocity",&
6587 hor_grid='v', z_grid='1')
6588 call register_restart_pair(cs%ubt_IC, cs%vbt_IC, vd(2), vd(3), .false., restart_cs, &
6589 conversion=us%L_T_to_m_s)
6590 endif
6591
6592 call register_restart_field(cs%dtbt, "DTBT", .false., restart_cs, &
6593 longname="Barotropic timestep", units="seconds", conversion=us%T_to_s)
6594
6595 ! Register streaming band-pass filters
6596 call get_param(param_file, mdl, "USE_FILTER", cs%use_filter, &
6597 "If true, use streaming band-pass filters to detect the "//&
6598 "instantaneous tidal signals in the simulation.", default=.false.)
6599 call get_param(param_file, mdl, "N_FILTERS", n_filters, &
6600 "Number of streaming band-pass filters to be used in the simulation.", &
6601 default=0, do_not_log=.not.cs%use_filter)
6602 if (n_filters<=0) cs%use_filter = .false.
6603 if (cs%use_filter) then
6604 call filt_register(n_filters, 'ubt', 'u', hi, cs%Filt_CS_u, restart_cs)
6605 call filt_register(n_filters, 'vbt', 'v', hi, cs%Filt_CS_v, restart_cs)
6606 endif
6607
6608end subroutine register_barotropic_restarts
6609
6610
6611!> \namespace mom_barotropic
6612!!
6613!! By Robert Hallberg, April 1994 - January 2007
6614!!
6615!! This program contains the subroutines that time steps the
6616!! linearized barotropic equations. btstep is used to actually
6617!! time step the barotropic equations, and contains most of the
6618!! substance of this module.
6619!!
6620!! btstep uses a forwards-backwards based scheme to time step
6621!! the barotropic equations, returning the layers' accelerations due
6622!! to the barotropic changes in the ocean state, the final free
6623!! surface height (or column mass), and the volume (or mass) fluxes
6624!! summed through the layers and averaged over the baroclinic time
6625!! step. As input, btstep takes the initial 3-D velocities, the
6626!! initial free surface height, the 3-D accelerations of the layers,
6627!! and the external forcing. Everything in btstep is cast in terms
6628!! of anomalies, so if everything is in balance, there is explicitly
6629!! no acceleration due to btstep.
6630!!
6631!! The spatial discretization of the continuity equation is second
6632!! order accurate. A flux conservative form is used to guarantee
6633!! global conservation of volume. The spatial discretization of the
6634!! momentum equation is second order accurate. The Coriolis force
6635!! is written in a form which does not contribute to the energy
6636!! tendency and which conserves linearized potential vorticity, f/D.
6637!! These terms are exactly removed from the baroclinic momentum
6638!! equations, so the linearization of vorticity advection will not
6639!! degrade the overall solution.
6640!!
6641!! btcalc calculates the fractional thickness of each layer at the
6642!! velocity points, for later use in calculating the barotropic
6643!! velocities and the averaged accelerations. Harmonic mean
6644!! thicknesses (i.e. 2*h_L*h_R/(h_L + h_R)) are used to avoid overly
6645!! strong weighting of overly thin layers. This may later be relaxed
6646!! to use thicknesses determined from the continuity equations.
6647!!
6648!! bt_mass_source determines the real mass sources for the
6649!! barotropic solver, along with the corrective pseudo-fluxes that
6650!! keep the barotropic and baroclinic estimates of the free surface
6651!! height close to each other. Given the layer thicknesses and the
6652!! free surface height that correspond to each other, it calculates
6653!! a corrective mass source that is added to the barotropic continuity*
6654!! equation, and optionally adjusts a slowly varying correction rate.
6655!! Newer algorithmic changes have deemphasized the need for this, but
6656!! it is still here to add net water sources to the barotropic solver.*
6657!!
6658!! barotropic_init allocates and initializes any barotropic arrays
6659!! that have not been read from a restart file, reads parameters from
6660!! the inputfile, and sets up diagnostic fields.
6661!!
6662!! barotropic_end deallocates anything allocated in barotropic_init
6663!! or register_barotropic_restarts.
6664!!
6665!! register_barotropic_restarts is used to indicate any fields that
6666!! are private to the barotropic solver that need to be included in
6667!! the restart files, and to ensure that they are read.
6668
6669end module mom_barotropic