MOM6
MOM_offline_main.F90
1 !> The routines here implement the offline tracer algorithm used in MOM6. These are called from step_offline
2 !! Some routines called here can be found in the MOM_offline_aux module.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use mpp_domains_mod, only : center, corner, north, east
8 use mom_ale, only : ale_cs, ale_main_offline, ale_offline_inputs, ale_offline_tracer_final
9 use mom_checksums, only : hchksum, uvchksum
10 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
11 use mom_cpu_clock, only : clock_component, clock_subcomponent
12 use mom_cpu_clock, only : clock_module_driver, clock_module, clock_routine
13 use mom_diabatic_aux, only : diabatic_aux_cs, set_pen_shortwave
14 use mom_diabatic_driver, only : diabatic_cs, extract_diabatic_member
15 use mom_diabatic_aux, only : tridiagts
16 use mom_diag_mediator, only : diag_ctrl, post_data, register_diag_field
17 use mom_domains, only : sum_across_pes, pass_var, pass_vector
18 use mom_error_handler, only : mom_error, mom_mesg, fatal, warning
19 use mom_error_handler, only : calltree_enter, calltree_leave
21 use mom_forcing_type, only : forcing
22 use mom_grid, only : ocean_grid_type
24 use mom_offline_aux, only : update_offline_from_arrays, update_offline_from_files
25 use mom_offline_aux, only : next_modulo_time, offline_add_diurnal_sw
26 use mom_offline_aux, only : update_h_horizontal_flux, update_h_vertical_flux, limit_mass_flux_3d
27 use mom_offline_aux, only : distribute_residual_uh_barotropic, distribute_residual_vh_barotropic
28 use mom_offline_aux, only : distribute_residual_uh_upwards, distribute_residual_vh_upwards
31 use mom_time_manager, only : time_type
32 use mom_tracer_advect, only : tracer_advect_cs, advect_tracer
33 use mom_tracer_diabatic, only : applytracerboundaryfluxesinout
34 use mom_tracer_flow_control, only : tracer_flow_control_cs, call_tracer_column_fns, call_tracer_stocks
35 use mom_tracer_registry, only : tracer_registry_type, mom_tracer_chksum, mom_tracer_chkinv
38 
39 implicit none ; private
40 
41 #include "MOM_memory.h"
42 #include "version_variable.h"
43 
44 !> The control structure for the offline transport module
45 type, public :: offline_transport_cs ; private
46 
47  ! Pointers to relevant fields from the main MOM control structure
48  type(ale_cs), pointer :: ale_csp => null()
49  !< A pointer to the ALE control structure
50  type(diabatic_cs), pointer :: diabatic_csp => null()
51  !< A pointer to the diabatic control structure
52  type(diag_ctrl), pointer :: diag => null()
53  !< Structure that regulates diagnostic output
54  type(ocean_obc_type), pointer :: obc => null()
55  !< A pointer to the open boundary condition control structure
56  type(tracer_advect_cs), pointer :: tracer_adv_csp => null()
57  !< A pointer to the tracer advection control structure
58  type(opacity_cs), pointer :: opacity_csp => null()
59  !< A pointer to the opacity control structure
60  type(tracer_flow_control_cs), pointer :: tracer_flow_csp => null()
61  !< A pointer to control structure that orchestrates the calling of tracer packages
62  type(tracer_registry_type), pointer :: tracer_reg => null()
63  !< A pointer to the tracer registry
64  type(thermo_var_ptrs), pointer :: tv => null()
65  !< A structure pointing to various thermodynamic variables
66  type(ocean_grid_type), pointer :: g => null()
67  !< Pointer to a structure containing metrics and related information
68  type(verticalgrid_type), pointer :: gv => null()
69  !< Pointer to structure containing information about the vertical grid
70  type(optics_type), pointer :: optics => null()
71  !< Pointer to the optical properties type
72  type(diabatic_aux_cs), pointer :: diabatic_aux_csp => null()
73  !< Pointer to the diabatic_aux control structure
74 
75  !> Variables related to reading in fields from online run
76  integer :: start_index !< Timelevel to start
77  integer :: iter_no !< Timelevel to start
78  integer :: numtime !< How many timelevels in the input fields
79  integer :: accumulated_time !< Length of time accumulated in the current offline interval
80  ! Index of each of the variables to be read in with separate indices for each variable if they
81  ! are set off from each other in time
82  integer :: ridx_sum = -1 !< Read index offset of the summed variables
83  integer :: ridx_snap = -1 !< Read index offset of the snapshot variables
84  integer :: nk_input !< Number of input levels in the input fields
85  character(len=200) :: offlinedir !< Directory where offline fields are stored
86  character(len=200) :: & ! Names of input files
87  surf_file, & !< Contains surface fields (2d arrays)
88  snap_file, & !< Snapshotted fields (layer thicknesses)
89  sum_file, & !< Fields which are accumulated over time
90  mean_file !< Fields averaged over time
91  character(len=20) :: redistribute_method !< 'barotropic' if evenly distributing extra flow
92  !! throughout entire watercolumn, 'upwards',
93  !! if trying to do it just in the layers above
94  !! 'both' if both methods are used
95  character(len=20) :: mld_var_name !< Name of the mixed layer depth variable to use
96  logical :: fields_are_offset !< True if the time-averaged fields and snapshot fields are
97  !! offset by one time level
98  logical :: x_before_y !< Which horizontal direction is advected first
99  logical :: print_adv_offline !< Prints out some updates each advection sub interation
100  logical :: skip_diffusion !< Skips horizontal diffusion of tracers
101  logical :: read_sw !< Read in averaged values for shortwave radiation
102  logical :: read_mld !< Check to see whether mixed layer depths should be read in
103  logical :: diurnal_sw !< Adds a synthetic diurnal cycle on shortwave radiation
104  logical :: debug !< If true, write verbose debugging messages
105  logical :: redistribute_barotropic !< Redistributes column-summed residual transports throughout
106  !! a column weighted by thickness
107  logical :: redistribute_upwards !< Redistributes remaining fluxes only in layers above
108  !! the current one based as the max allowable transport
109  !! in that cell
110  logical :: read_all_ts_uvh !< If true, then all timelevels of temperature, salinity, mass transports, and
111  !! Layer thicknesses are read during initialization
112  !! Variables controlling some of the numerical considerations of offline transport
113  integer :: num_off_iter !< Number of advection iterations per offline step
114  integer :: num_vert_iter !< Number of vertical iterations per offline step
115  integer :: off_ale_mod !< Sets how frequently the ALE step is done during the advection
116  real :: dt_offline !< Timestep used for offline tracers [s]
117  real :: dt_offline_vertical !< Timestep used for calls to tracer vertical physics [s]
118  real :: evap_cfl_limit !< Copied from diabatic_CS controlling how tracers follow freshwater fluxes
119  real :: minimum_forcing_depth !< Copied from diabatic_CS controlling how tracers follow freshwater fluxes
120  real :: kd_max !< Runtime parameter specifying the maximum value of vertical diffusivity
121  real :: min_residual !< The minimum amount of total mass flux before exiting the main advection routine
122  !>@{ Diagnostic manager IDs for some fields that may be of interest when doing offline transport
123  integer :: &
124  id_uhr = -1, &
125  id_vhr = -1, &
126  id_ear = -1, &
127  id_ebr = -1, &
128  id_hr = -1, &
129  id_hdiff = -1, &
130  id_uhr_redist = -1, &
131  id_vhr_redist = -1, &
132  id_uhr_end = -1, &
133  id_vhr_end = -1, &
134  id_eta_pre_distribute = -1, &
135  id_eta_post_distribute = -1, &
136  id_h_redist = -1, &
137  id_eta_diff_end = -1
138 
139  ! Diagnostic IDs for the regridded/remapped input fields
140  integer :: &
141  id_uhtr_regrid = -1, &
142  id_vhtr_regrid = -1, &
143  id_temp_regrid = -1, &
144  id_salt_regrid = -1, &
145  id_h_regrid = -1
146  !!@}
147 
148  ! IDs for timings of various offline components
149  integer :: id_clock_read_fields = -1 !< A CPU time clock
150  integer :: id_clock_offline_diabatic = -1 !< A CPU time clock
151  integer :: id_clock_offline_adv = -1 !< A CPU time clock
152  integer :: id_clock_redistribute = -1 !< A CPU time clock
153 
154  !> Zonal transport that may need to be stored between calls to step_MOM
155  real, allocatable, dimension(:,:,:) :: uhtr
156  !> Meridional transport that may need to be stored between calls to step_MOM
157  real, allocatable, dimension(:,:,:) :: vhtr
158 
159  ! Fields at T-point
160  real, allocatable, dimension(:,:,:) :: eatr
161  !< Amount of fluid entrained from the layer above within
162  !! one time step [H ~> m or kg m-2]
163  real, allocatable, dimension(:,:,:) :: ebtr
164  !< Amount of fluid entrained from the layer below within
165  !! one time step [H ~> m or kg m-2]
166  ! Fields at T-points on interfaces
167  real, allocatable, dimension(:,:,:) :: kd !< Vertical diffusivity
168  real, allocatable, dimension(:,:,:) :: h_end !< Thicknesses at the end of offline timestep
169 
170  real, allocatable, dimension(:,:) :: netmassin !< Freshwater fluxes into the ocean
171  real, allocatable, dimension(:,:) :: netmassout !< Freshwater fluxes out of the ocean
172  real, allocatable, dimension(:,:) :: mld !< Mixed layer depths at thickness points [H ~> m or kg m-2].
173 
174  ! Allocatable arrays to read in entire fields during initialization
175  real, allocatable, dimension(:,:,:,:) :: uhtr_all !< Entire field of zonal transport
176  real, allocatable, dimension(:,:,:,:) :: vhtr_all !< Entire field of mericional transport
177  real, allocatable, dimension(:,:,:,:) :: hend_all !< Entire field of layer thicknesses
178  real, allocatable, dimension(:,:,:,:) :: temp_all !< Entire field of temperatures
179  real, allocatable, dimension(:,:,:,:) :: salt_all !< Entire field of salinities
180 
181 end type offline_transport_cs
182 
183 public offline_advection_ale
184 public offline_redistribute_residual
185 public offline_diabatic_ale
186 public offline_fw_fluxes_into_ocean
187 public offline_fw_fluxes_out_ocean
188 public offline_advection_layer
189 public register_diags_offline_transport
190 public update_offline_fields
191 public insert_offline_main
192 public extract_offline_main
193 public post_offline_convergence_diags
194 public offline_transport_init
195 public offline_transport_end
196 
197 contains
198 
199 !> 3D advection is done by doing flux-limited nonlinear horizontal advection interspersed with an ALE
200 !! regridding/remapping step. The loop in this routine is exited if remaining residual transports are below
201 !! a runtime-specified value or a maximum number of iterations is reached.
202 subroutine offline_advection_ale(fluxes, Time_start, time_interval, CS, id_clock_ale, h_pre, uhtr, vhtr, converged)
203  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
204  type(time_type), intent(in) :: time_start !< starting time of a segment, as a time type
205  real, intent(in) :: time_interval !< time interval
206  type(offline_transport_cs), pointer :: cs !< control structure for offline module
207  integer, intent(in) :: id_clock_ale !< Clock for ALE routines
208  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
209  intent(inout) :: h_pre !< layer thicknesses before advection
210  !! [H ~> m or kg m-2]
211  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
212  intent(inout) :: uhtr !< Zonal mass transport [H m2 ~> m3 or kg]
213  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), &
214  intent(inout) :: vhtr !< Meridional mass transport [H m2 ~> m3 or kg]
215  logical, intent( out) :: converged !< True if the iterations have converged
216 
217  ! Local pointers
218  type(ocean_grid_type), pointer :: g => null() ! Pointer to a structure containing
219  ! metrics and related information
220  type(verticalgrid_type), pointer :: gv => null() ! Pointer to structure containing information
221  ! about the vertical grid
222  ! Work arrays for mass transports
223  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhtr_sub
224  ! Meridional mass transports
225  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhtr_sub
226 
227  real :: prev_tot_residual, tot_residual ! Used to keep track of how close to convergence we are
228 
229  ! Variables used to keep track of layer thicknesses at various points in the code
230  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
231  h_new, &
232  h_vol
233  ! Fields for eta_diff diagnostic
234  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_pre, eta_end
235  integer :: niter, iter
236  real :: inum_iter
237  character(len=256) :: mesg ! The text of an error message
238  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz
239  integer :: isv, iev, jsv, jev ! The valid range of the indices.
240  integer :: isdb, iedb, jsdb, jedb
241  logical :: z_first, x_before_y
242  real :: evap_cfl_limit, minimum_forcing_depth, dt_iter, dt_offline
243 
244  integer :: nstocks
245  real :: stock_values(max_fields_)
246  character(len=20) :: debug_msg
247  call cpu_clock_begin(cs%id_clock_offline_adv)
248 
249  ! Grid-related pointer assignments
250  g => cs%G
251  gv => cs%GV
252 
253  x_before_y = cs%x_before_y
254 
255  ! Initialize some shorthand variables from other structures
256  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
257  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
258  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
259 
260  dt_offline = cs%dt_offline
261  evap_cfl_limit = cs%evap_CFL_limit
262  minimum_forcing_depth = cs%minimum_forcing_depth
263 
264  niter = cs%num_off_iter
265  inum_iter = 1./real(niter)
266  dt_iter = dt_offline*inum_iter
267 
268  ! Initialize working arrays
269  h_new(:,:,:) = 0.0
270  h_vol(:,:,:) = 0.0
271  uhtr_sub(:,:,:) = 0.0
272  vhtr_sub(:,:,:) = 0.0
273 
274  ! converged should only be true if there are no remaining mass fluxes
275  converged = .false.
276 
277  ! Tracers are transported using the stored mass fluxes. Where possible, operators are Strang-split around
278  ! the call to
279  ! 1) Using the layer thicknesses and tracer concentrations from the previous timestep,
280  ! half of the accumulated vertical mixing (eatr and ebtr) is applied in the call to tracer_column_fns.
281  ! For tracers whose source/sink terms need dt, this value is set to 1/2 dt_offline
282  ! 2) Half of the accumulated surface freshwater fluxes are applied
283  !! START ITERATION
284  ! 3) Accumulated mass fluxes are used to do horizontal transport. The number of iterations used in
285  ! advect_tracer is limited to 2 (e.g x->y->x->y). The remaining mass fluxes are stored for later use
286  ! and resulting layer thicknesses fed into the next step
287  ! 4) Tracers and the h-grid are regridded and remapped in a call to ALE. This allows for layers which might
288  ! 'vanish' because of horizontal mass transport to be 'reinflated'
289  ! 5) Check that transport is done if the remaining mass fluxes equals 0 or if the max number of iterations
290  ! has been reached
291  !! END ITERATION
292  ! 6) Repeat steps 1 and 2
293  ! 7) Force a remapping to the stored layer thicknesses that correspond to the snapshot of the online model
294  ! 8) Reset T/S and h to their stored snapshotted values to prevent model drift
295 
296  ! Copy over the horizontal mass fluxes from the total mass fluxes
297  do k=1,nz ; do j=jsd,jed ; do i=isdb,iedb
298  uhtr_sub(i,j,k) = uhtr(i,j,k)
299  enddo ; enddo ; enddo
300  do k=1,nz ; do j=jsdb,jedb ; do i=isd,ied
301  vhtr_sub(i,j,k) = vhtr(i,j,k)
302  enddo ; enddo ; enddo
303  do k=1,nz ; do j=js,je ; do i=is,ie
304  h_new(i,j,k) = h_pre(i,j,k)
305  enddo ; enddo ; enddo
306 
307  if (cs%debug) then
308  call hchksum(h_pre,"h_pre before transport",g%HI)
309  call uvchksum("[uv]htr_sub before transport", uhtr_sub, vhtr_sub, g%HI)
310  endif
311  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
312  if (cs%print_adv_offline) then
313  write(mesg,'(A,ES24.16)') "Main advection starting transport: ", tot_residual
314  call mom_mesg(mesg)
315  endif
316 
317  ! This loop does essentially a flux-limited, nonlinear advection scheme until all mass fluxes
318  ! are used. ALE is done after the horizontal advection.
319  do iter=1,cs%num_off_iter
320 
321  do k=1,nz ; do j=js,je ; do i=is,ie
322  h_vol(i,j,k) = h_new(i,j,k)*g%areaT(i,j)
323  h_pre(i,j,k) = h_new(i,j,k)
324  enddo ; enddo ; enddo
325 
326  if (cs%debug) then
327  call hchksum(h_vol,"h_vol before advect",g%HI)
328  call uvchksum("[uv]htr_sub before advect", uhtr_sub, vhtr_sub, g%HI)
329  write(debug_msg, '(A,I4.4)') 'Before advect ', iter
330  call mom_tracer_chkinv(debug_msg, g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
331  endif
332 
333  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, cs%dt_offline, g, gv, &
334  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=1, &
335  uhr_out=uhtr, vhr_out=vhtr, h_out=h_new, x_first_in=x_before_y)
336 
337  ! Switch the direction every iteration
338  x_before_y = .not. x_before_y
339 
340  ! Update the new layer thicknesses after one round of advection has happened
341  do k=1,nz ; do j=js,je ; do i=is,ie
342  h_new(i,j,k) = h_new(i,j,k)/g%areaT(i,j)
343  enddo ; enddo ; enddo
344 
345  if (modulo(iter,cs%off_ale_mod)==0) then
346  ! Do ALE remapping/regridding to allow for more advection to occur in the next iteration
347  call pass_var(h_new,g%Domain)
348  if (cs%debug) then
349  call hchksum(h_new,"h_new before ALE",g%HI)
350  write(debug_msg, '(A,I4.4)') 'Before ALE ', iter
351  call mom_tracer_chkinv(debug_msg, g, h_new, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
352  endif
353  call cpu_clock_begin(id_clock_ale)
354  call ale_main_offline(g, gv, h_new, cs%tv, cs%tracer_Reg, cs%ALE_CSp, cs%dt_offline)
355  call cpu_clock_end(id_clock_ale)
356 
357  if (cs%debug) then
358  call hchksum(h_new,"h_new after ALE",g%HI)
359  write(debug_msg, '(A,I4.4)') 'After ALE ', iter
360  call mom_tracer_chkinv(debug_msg, g, h_new, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
361  endif
362  endif
363 
364  do k=1,nz; do j=js,je ; do i=is,ie
365  uhtr_sub(i,j,k) = uhtr(i,j,k)
366  vhtr_sub(i,j,k) = vhtr(i,j,k)
367  enddo ; enddo ; enddo
368  call pass_var(h_new, g%Domain)
369  call pass_vector(uhtr_sub,vhtr_sub,g%Domain)
370 
371  ! Check for whether we've used up all the advection, or if we need to move on because
372  ! advection has stalled
373  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
374  if (cs%print_adv_offline) then
375  write(mesg,'(A,ES24.16)') "Main advection remaining transport: ", tot_residual
376  call mom_mesg(mesg)
377  endif
378  ! If all the mass transports have been used u, then quit
379  if (tot_residual == 0.0) then
380  write(mesg,*) "Converged after iteration ", iter
381  call mom_mesg(mesg)
382  converged = .true.
383  exit
384  endif
385  ! If advection has stalled or the remaining residual is less than a specified amount, quit
386  if ( (tot_residual == prev_tot_residual) .or. (tot_residual<cs%min_residual) ) then
387  converged = .false.
388  exit
389  endif
390 
391  prev_tot_residual = tot_residual
392 
393  enddo
394 
395  ! Make sure that uhtr and vhtr halos are updated
396  h_pre(:,:,:) = h_new(:,:,:)
397  call pass_vector(uhtr,vhtr,g%Domain)
398 
399  if (cs%debug) then
400  call hchksum(h_pre,"h after offline_advection_ale",g%HI)
401  call uvchksum("[uv]htr after offline_advection_ale", uhtr, vhtr, g%HI)
402  call mom_tracer_chkinv("After offline_advection_ale", g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
403  endif
404 
405  call cpu_clock_end(cs%id_clock_offline_adv)
406 
407 end subroutine offline_advection_ale
408 
409 !> In the case where the main advection routine did not converge, something needs to be done with the remaining
410 !! transport. Two different ways are offered, 'barotropic' means that the residual is distributed equally
411 !! throughout the water column. 'upwards' attempts to redistribute the transport in the layers above and will
412 !! eventually work down the entire water column
413 subroutine offline_redistribute_residual(CS, h_pre, uhtr, vhtr, converged)
414  type(offline_transport_cs), pointer :: cs !< control structure from initialize_MOM
415  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
416  intent(inout) :: h_pre !< layer thicknesses before advection
417  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
418  intent(inout) :: uhtr !< Zonal mass transport
419  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), &
420  intent(inout) :: vhtr !< Meridional mass transport
421  logical, intent(in ) :: converged !< True if the iterations have converged
422 
423  type(ocean_grid_type), pointer :: g => null() ! Pointer to a structure containing
424  ! metrics and related information
425  type(verticalgrid_type), pointer :: gv => null() ! Pointer to structure containing information
426  ! about the vertical grid
427  logical :: x_before_y
428  ! Variables used to keep track of layer thicknesses at various points in the code
429  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
430  h_new, &
431  h_vol
432 
433  ! Used to calculate the eta diagnostics
434  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_work
435  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhr !< Zonal mass transport
436  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhr !< Meridional mass transport
437 
438  character(len=256) :: mesg ! The text of an error message
439  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz, iter
440  real :: prev_tot_residual, tot_residual, stock_values(max_fields_)
441  integer :: nstocks
442 
443  ! Assign grid pointers
444  g => cs%G
445  gv => cs%GV
446 
447  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
448  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
449 
450  x_before_y = cs%x_before_y
451 
452  if (cs%id_eta_pre_distribute>0) then
453  eta_work(:,:) = 0.0
454  do k=1,nz ; do j=js,je ; do i=is,ie
455  if (h_pre(i,j,k)>gv%Angstrom_H) then
456  eta_work(i,j) = eta_work(i,j) + h_pre(i,j,k)
457  endif
458  enddo ; enddo ; enddo
459  call post_data(cs%id_eta_pre_distribute,eta_work,cs%diag)
460  endif
461 
462  ! These are used to find out how much will be redistributed in this routine
463  if (cs%id_h_redist>0) call post_data(cs%id_h_redist, h_pre, cs%diag)
464  if (cs%id_uhr_redist>0) call post_data(cs%id_uhr_redist, uhtr, cs%diag)
465  if (cs%id_vhr_redist>0) call post_data(cs%id_vhr_redist, vhtr, cs%diag)
466 
467  if (converged) return
468 
469  if (cs%debug) then
470  call mom_tracer_chkinv("Before redistribute ", g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
471  endif
472 
473  call cpu_clock_begin(cs%id_clock_redistribute)
474 
475  if (cs%redistribute_upwards .or. cs%redistribute_barotropic) then
476  do iter = 1, cs%num_off_iter
477 
478  ! Perform upwards redistribution
479  if (cs%redistribute_upwards) then
480 
481  ! Calculate the layer volumes at beginning of redistribute
482  do k=1,nz ; do j=js,je ; do i=is,ie
483  h_vol(i,j,k) = h_pre(i,j,k)*g%areaT(i,j)
484  enddo ; enddo ; enddo
485  call pass_var(h_vol,g%Domain)
486  call pass_vector(uhtr,vhtr,g%Domain)
487 
488  ! Store volumes for advect_tracer
489  h_pre(:,:,:) = h_vol(:,:,:)
490 
491  if (cs%debug) then
492  call mom_tracer_chksum("Before upwards redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
493  call uvchksum("[uv]tr before upwards redistribute", uhtr, vhtr, g%HI)
494  endif
495 
496  if (x_before_y) then
497  call distribute_residual_uh_upwards(g, gv, h_vol, uhtr)
498  call distribute_residual_vh_upwards(g, gv, h_vol, vhtr)
499  else
500  call distribute_residual_vh_upwards(g, gv, h_vol, vhtr)
501  call distribute_residual_uh_upwards(g, gv, h_vol, uhtr)
502  endif
503 
504  call advect_tracer(h_pre, uhtr, vhtr, cs%OBC, cs%dt_offline, g, gv, &
505  cs%tracer_adv_CSp, cs%tracer_Reg, h_prev_opt = h_pre, max_iter_in=1, &
506  h_out=h_new, uhr_out=uhr, vhr_out=vhr, x_first_in=x_before_y)
507 
508  if (cs%debug) then
509  call mom_tracer_chksum("After upwards redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
510  endif
511 
512  ! Convert h_new back to layer thickness for ALE remapping
513  do k=1,nz ; do j=js,je ; do i=is,ie
514  uhtr(i,j,k) = uhr(i,j,k)
515  vhtr(i,j,k) = vhr(i,j,k)
516  h_vol(i,j,k) = h_new(i,j,k)
517  h_new(i,j,k) = h_new(i,j,k)/g%areaT(i,j)
518  h_pre(i,j,k) = h_new(i,j,k)
519  enddo ; enddo ; enddo
520 
521  endif ! redistribute upwards
522 
523  ! Perform barotropic redistribution
524  if (cs%redistribute_barotropic) then
525 
526  ! Calculate the layer volumes at beginning of redistribute
527  do k=1,nz ; do j=js,je ; do i=is,ie
528  h_vol(i,j,k) = h_pre(i,j,k)*g%areaT(i,j)
529  enddo ; enddo ; enddo
530  call pass_var(h_vol,g%Domain)
531  call pass_vector(uhtr,vhtr,g%Domain)
532 
533  ! Copy h_vol to h_pre for advect_tracer routine
534  h_pre(:,:,:) = h_vol(:,:,:)
535 
536  if (cs%debug) then
537  call mom_tracer_chksum("Before barotropic redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
538  call uvchksum("[uv]tr before upwards redistribute", uhtr, vhtr, g%HI)
539  endif
540 
541  if (x_before_y) then
542  call distribute_residual_uh_barotropic(g, gv, h_vol, uhtr)
543  call distribute_residual_vh_barotropic(g, gv, h_vol, vhtr)
544  else
545  call distribute_residual_vh_barotropic(g, gv, h_vol, vhtr)
546  call distribute_residual_uh_barotropic(g, gv, h_vol, uhtr)
547  endif
548 
549  call advect_tracer(h_pre, uhtr, vhtr, cs%OBC, cs%dt_offline, g, gv, &
550  cs%tracer_adv_CSp, cs%tracer_Reg, h_prev_opt = h_pre, max_iter_in=1, &
551  h_out=h_new, uhr_out=uhr, vhr_out=vhr, x_first_in=x_before_y)
552 
553  if (cs%debug) then
554  call mom_tracer_chksum("After barotropic redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
555  endif
556 
557  ! Convert h_new back to layer thickness for ALE remapping
558  do k=1,nz ; do j=js,je ; do i=is,ie
559  uhtr(i,j,k) = uhr(i,j,k)
560  vhtr(i,j,k) = vhr(i,j,k)
561  h_vol(i,j,k) = h_new(i,j,k)
562  h_new(i,j,k) = h_new(i,j,k)/g%areaT(i,j)
563  h_pre(i,j,k) = h_new(i,j,k)
564  enddo ; enddo ; enddo
565 
566  endif ! redistribute barotropic
567 
568  ! Check to see if all transport has been exhausted
569  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
570  if (cs%print_adv_offline) then
571  write(mesg,'(A,ES24.16)') "Residual advection remaining transport: ", tot_residual
572  call mom_mesg(mesg)
573  endif
574  ! If the remaining residual is 0, then this return is done
575  if (tot_residual==0.0 ) then
576  exit
577  endif
578 
579  if ( (tot_residual == prev_tot_residual) .or. (tot_residual<cs%min_residual) ) exit
580  prev_tot_residual = tot_residual
581 
582  enddo ! Redistribution iteration
583  endif ! If one of the redistribution routines is requested
584 
585  if (cs%id_eta_post_distribute>0) then
586  eta_work(:,:) = 0.0
587  do k=1,nz ; do j=js,je ; do i=is,ie
588  if (h_pre(i,j,k)>gv%Angstrom_H) then
589  eta_work(i,j) = eta_work(i,j) + h_pre(i,j,k)
590  endif
591  enddo ; enddo ; enddo
592  call post_data(cs%id_eta_post_distribute,eta_work,cs%diag)
593  endif
594 
595  if (cs%id_uhr>0) call post_data(cs%id_uhr,uhtr,cs%diag)
596  if (cs%id_vhr>0) call post_data(cs%id_vhr,vhtr,cs%diag)
597 
598  if (cs%debug) then
599  call hchksum(h_pre,"h_pre after redistribute",g%HI)
600  call uvchksum("uhtr after redistribute", uhtr, vhtr, g%HI)
601  call mom_tracer_chkinv("after redistribute ", g, h_new, cs%tracer_Reg%Tr, cs%tracer_Reg%ntr)
602  endif
603 
604  call cpu_clock_end(cs%id_clock_redistribute)
605 
606 end subroutine offline_redistribute_residual
607 
608 !> Sums any non-negligible remaining transport to check for advection convergence
609 real function remaining_transport_sum(CS, uhtr, vhtr)
610  type(offline_transport_cs), pointer :: cs !< control structure for offline module
611  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(in ) :: uhtr !< Zonal mass transport
612  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(in ) :: vhtr !< Meridional mass transport
613 
614  ! Local variables
615  integer :: i, j, k
616  integer :: is, ie, js, je, nz
617  real :: h_min !< A layer thickness below roundoff from GV type
618  real :: uh_neglect !< A small value of zonal transport that effectively is below roundoff error
619  real :: vh_neglect !< A small value of meridional transport that effectively is below roundoff error
620 
621  nz = cs%GV%ke
622  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec
623 
624  h_min = cs%GV%H_subroundoff
625 
626  remaining_transport_sum = 0.
627  do k=1,nz; do j=js,je ; do i=is,ie
628  uh_neglect = h_min*min(cs%G%areaT(i,j),cs%G%areaT(i+1,j))
629  vh_neglect = h_min*min(cs%G%areaT(i,j),cs%G%areaT(i,j+1))
630  if (abs(uhtr(i,j,k))>uh_neglect) then
631  remaining_transport_sum = remaining_transport_sum + abs(uhtr(i,j,k))
632  endif
633  if (abs(vhtr(i,j,k))>vh_neglect) then
634  remaining_transport_sum = remaining_transport_sum + abs(vhtr(i,j,k))
635  endif
636  enddo ; enddo ; enddo
637  call sum_across_pes(remaining_transport_sum)
638 
639 end function remaining_transport_sum
640 
641 !> The vertical/diabatic driver for offline tracers. First the eatr/ebtr associated with the interpolated
642 !! vertical diffusivities are calculated and then any tracer column functions are done which can include
643 !! vertical diffuvities and source/sink terms.
644 subroutine offline_diabatic_ale(fluxes, Time_start, Time_end, CS, h_pre, eatr, ebtr)
645 
646  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
647  type(time_type), intent(in) :: time_start !< starting time of a segment, as a time type
648  type(time_type), intent(in) :: time_end !< time interval
649  type(offline_transport_cs), pointer :: cs !< control structure from initialize_MOM
650  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
651  intent(inout) :: h_pre !< layer thicknesses before advection [H ~> m or kg m-2]
652  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
653  intent(inout) :: eatr !< Entrainment from layer above [H ~> m or kg m-2]
654  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
655  intent(inout) :: ebtr !< Entrainment from layer below [H ~> m or kg m-2]
656 
657  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: sw, sw_vis, sw_nir !< Save old value of shortwave radiation
658  real :: hval
659  integer :: i,j,k
660  integer :: is, ie, js, je, nz
661  integer :: k_nonzero
662  real :: stock_values(max_fields_)
663  real :: kd_bot
664  integer :: nstocks
665  nz = cs%GV%ke
666  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec
667 
668  call cpu_clock_begin(cs%id_clock_offline_diabatic)
669 
670  call mom_mesg("Applying tracer source, sinks, and vertical mixing")
671 
672  if (cs%debug) then
673  call hchksum(h_pre,"h_pre before offline_diabatic_ale",cs%G%HI)
674  call hchksum(eatr,"eatr before offline_diabatic_ale",cs%G%HI)
675  call hchksum(ebtr,"ebtr before offline_diabatic_ale",cs%G%HI)
676  call mom_tracer_chkinv("Before offline_diabatic_ale", cs%G, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
677  endif
678 
679  eatr(:,:,:) = 0.
680  ebtr(:,:,:) = 0.
681  ! Calculate eatr and ebtr if vertical diffusivity is read
682  ! Because the saved remapped diagnostics from the online run assume a zero minimum thickness
683  ! but ALE may have a minimum thickness. Flood the diffusivities for all layers with the value
684  ! of Kd closest to the bottom which is non-zero
685  do j=js,je ; do i=is,ie
686  k_nonzero = nz+1
687  ! Find the nonzero bottom Kd
688  do k=nz+1,1,-1
689  if (cs%Kd(i,j,k)>0.) then
690  kd_bot = cs%Kd(i,j,k)
691  k_nonzero = k
692  exit
693  endif
694  enddo
695  ! Flood the bottom interfaces
696  do k=k_nonzero,nz+1
697  cs%Kd(i,j,k) = kd_bot
698  enddo
699  enddo ; enddo
700 
701  do j=js,je ; do i=is,ie
702  eatr(i,j,1) = 0.
703  enddo ; enddo
704  do k=2,nz ; do j=js,je ; do i=is,ie
705  hval=1.0/(cs%GV%H_subroundoff + 0.5*(h_pre(i,j,k-1) + h_pre(i,j,k)))
706  eatr(i,j,k) = (cs%GV%m_to_H**2) * cs%dt_offline_vertical * hval * cs%Kd(i,j,k)
707  ebtr(i,j,k-1) = eatr(i,j,k)
708  enddo ; enddo ; enddo
709  do j=js,je ; do i=is,ie
710  ebtr(i,j,nz) = 0.
711  enddo ; enddo
712 
713  ! Add diurnal cycle for shortwave radiation (only used if run in ocean-only mode)
714  if (cs%diurnal_SW .and. cs%read_sw) then
715  sw(:,:) = fluxes%sw
716  sw_vis(:,:) = fluxes%sw_vis_dir
717  sw_nir(:,:) = fluxes%sw_nir_dir
718  call offline_add_diurnal_sw(fluxes, cs%G, time_start, time_end)
719  endif
720 
721  if (associated(cs%optics)) &
722  call set_pen_shortwave(cs%optics, fluxes, cs%G, cs%GV, cs%diabatic_aux_CSp, cs%opacity_CSp, cs%tracer_flow_CSp)
723 
724  ! Note that tracerBoundaryFluxesInOut within this subroutine should NOT be called
725  ! as the freshwater fluxes have already been accounted for
726  call call_tracer_column_fns(h_pre, h_pre, eatr, ebtr, fluxes, cs%MLD, cs%dt_offline_vertical, cs%G, cs%GV, &
727  cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
728 
729  if (cs%diurnal_SW .and. cs%read_sw) then
730  fluxes%sw(:,:) = sw
731  fluxes%sw_vis_dir(:,:) = sw_vis
732  fluxes%sw_nir_dir(:,:) = sw_nir
733  endif
734 
735  if (cs%debug) then
736  call hchksum(h_pre,"h_pre after offline_diabatic_ale",cs%G%HI)
737  call hchksum(eatr,"eatr after offline_diabatic_ale",cs%G%HI)
738  call hchksum(ebtr,"ebtr after offline_diabatic_ale",cs%G%HI)
739  call mom_tracer_chkinv("After offline_diabatic_ale", cs%G, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
740  endif
741 
742  call cpu_clock_end(cs%id_clock_offline_diabatic)
743 
744 end subroutine offline_diabatic_ale
745 
746 !> Apply positive freshwater fluxes (into the ocean) and update netMassOut with only the negative
747 !! (out of the ocean) fluxes
748 subroutine offline_fw_fluxes_into_ocean(G, GV, CS, fluxes, h, in_flux_optional)
749  type(offline_transport_cs), intent(inout) :: cs !< Offline control structure
750  type(ocean_grid_type), intent(in) :: g !< Grid structure
751  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
752  type(forcing), intent(inout) :: fluxes !< Surface fluxes container
753  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
754  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
755  real, dimension(SZI_(G),SZJ_(G)), &
756  optional, intent(in) :: in_flux_optional !< The total time-integrated amount
757  !! of tracer that leaves with freshwater
758 
759  integer :: i, j, m
760  real, dimension(SZI_(G),SZJ_(G)) :: negative_fw !< store all negative fluxes
761  logical :: update_h !< Flag for whether h should be updated
762 
763  if ( present(in_flux_optional) ) &
764  call mom_error(warning, "Positive freshwater fluxes with non-zero tracer concentration not supported yet")
765 
766  ! Set all fluxes to 0
767  negative_fw(:,:) = 0.
768 
769  ! Sort fluxes into positive and negative
770  do j=g%jsc,g%jec ; do i=g%isc,g%iec
771  if (fluxes%netMassOut(i,j)<0.0) then
772  negative_fw(i,j) = fluxes%netMassOut(i,j)
773  fluxes%netMassOut(i,j) = 0.
774  endif
775  enddo ; enddo
776 
777  if (cs%debug) then
778  call hchksum(h,"h before fluxes into ocean",g%HI)
779  call mom_tracer_chkinv("Before fluxes into ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
780  endif
781  do m = 1,cs%tracer_reg%ntr
782  ! Layer thicknesses should only be updated after the last tracer is finished
783  update_h = ( m == cs%tracer_reg%ntr )
784  call applytracerboundaryfluxesinout(g, gv, cs%tracer_reg%tr(m)%t, cs%dt_offline, fluxes, h, &
785  cs%evap_CFL_limit, cs%minimum_forcing_depth, update_h_opt = update_h)
786  enddo
787  if (cs%debug) then
788  call hchksum(h,"h after fluxes into ocean",g%HI)
789  call mom_tracer_chkinv("After fluxes into ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
790  endif
791 
792  ! Now that fluxes into the ocean are done, save the negative fluxes for later
793  fluxes%netMassOut(:,:) = negative_fw(:,:)
794 
795 end subroutine offline_fw_fluxes_into_ocean
796 
797 !> Apply negative freshwater fluxes (out of the ocean)
798 subroutine offline_fw_fluxes_out_ocean(G, GV, CS, fluxes, h, out_flux_optional)
799  type(offline_transport_cs), intent(inout) :: cs !< Offline control structure
800  type(ocean_grid_type), intent(in) :: g !< Grid structure
801  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
802  type(forcing), intent(inout) :: fluxes !< Surface fluxes container
803  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
804  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
805  real, dimension(SZI_(G),SZJ_(G)), &
806  optional, intent(in) :: out_flux_optional !< The total time-integrated amount
807  !! of tracer that leaves with freshwater
808 
809  integer :: m
810  logical :: update_h !< Flag for whether h should be updated
811 
812  if ( present(out_flux_optional) ) &
813  call mom_error(warning, "Negative freshwater fluxes with non-zero tracer concentration not supported yet")
814 
815  if (cs%debug) then
816  call hchksum(h,"h before fluxes out of ocean",g%HI)
817  call mom_tracer_chkinv("Before fluxes out of ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
818  endif
819  do m = 1, cs%tracer_reg%ntr
820  ! Layer thicknesses should only be updated after the last tracer is finished
821  update_h = ( m == cs%tracer_reg%ntr )
822  call applytracerboundaryfluxesinout(g, gv, cs%tracer_reg%tr(m)%t, cs%dt_offline, fluxes, h, &
823  cs%evap_CFL_limit, cs%minimum_forcing_depth, update_h_opt = update_h)
824  enddo
825  if (cs%debug) then
826  call hchksum(h,"h after fluxes out of ocean",g%HI)
827  call mom_tracer_chkinv("Before fluxes out of ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
828  endif
829 
830 end subroutine offline_fw_fluxes_out_ocean
831 
832 !> When in layer mode, 3D horizontal advection using stored mass fluxes must be used. Horizontal advection is
833 !! done via tracer_advect, whereas the vertical component is actually handled by vertdiff in tracer_column_fns
834 subroutine offline_advection_layer(fluxes, Time_start, time_interval, CS, h_pre, eatr, ebtr, uhtr, vhtr)
835  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
836  type(time_type), intent(in) :: time_start !< starting time of a segment, as a time type
837  real, intent(in) :: time_interval !< Offline transport time interval
838  type(offline_transport_cs), pointer :: cs !< Control structure for offline module
839  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_pre !< layer thicknesses before advection
840  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: eatr !< Entrainment from layer above
841  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: ebtr !< Entrainment from layer below
842  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: uhtr !< Zonal mass transport
843  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(inout) :: vhtr !< Meridional mass transport
844  ! Local pointers
845  type(ocean_grid_type), pointer :: g => null() ! Pointer to a structure containing
846  ! metrics and related information
847  type(verticalgrid_type), pointer :: gv => null() ! Pointer to structure containing information
848  ! about the vertical grid
849  ! Remaining zonal mass transports
850  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhtr_sub
851  ! Remaining meridional mass transports
852  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhtr_sub
853 
854  real :: sum_abs_fluxes, sum_u, sum_v ! Used to keep track of how close to convergence we are
855  real :: dt_offline
856 
857  ! Local variables
858  ! Vertical diffusion related variables
859  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
860  eatr_sub, &
861  ebtr_sub
862  ! Variables used to keep track of layer thicknesses at various points in the code
863  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
864  h_new, &
865  h_vol
866  ! Work arrays for temperature and salinity
867  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
868  temp_old, salt_old, &
869  temp_mean, salt_mean, &
870  zero_3dh !
871  integer :: niter, iter
872  real :: inum_iter, dt_iter
873  logical :: converged
874  character(len=160) :: mesg ! The text of an error message
875  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz
876  integer :: isv, iev, jsv, jev ! The valid range of the indices.
877  integer :: isdb, iedb, jsdb, jedb
878  logical :: z_first, x_before_y
879 
880  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
881  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
882  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
883 
884  do iter=1,cs%num_off_iter
885 
886  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
887  eatr_sub(i,j,k) = eatr(i,j,k)
888  ebtr_sub(i,j,k) = ebtr(i,j,k)
889  enddo ; enddo ; enddo
890 
891  do k = 1, nz ; do j=js-1,je+1 ; do i=is-2,ie+1
892  uhtr_sub(i,j,k) = uhtr(i,j,k)
893  enddo ; enddo ; enddo
894 
895  do k = 1, nz ; do j=js-2,je+1 ; do i=is-1,ie+1
896  vhtr_sub(i,j,k) = vhtr(i,j,k)
897  enddo ; enddo ; enddo
898 
899 
900  ! Calculate 3d mass transports to be used in this iteration
901  call limit_mass_flux_3d(g, gv, uhtr_sub, vhtr_sub, eatr_sub, ebtr_sub, h_pre)
902 
903  if (z_first) then
904  ! First do vertical advection
905  call update_h_vertical_flux(g, gv, eatr_sub, ebtr_sub, h_pre, h_new)
906  call call_tracer_column_fns(h_pre, h_new, eatr_sub, ebtr_sub, &
907  fluxes, cs%mld, dt_iter, g, gv, cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
908  ! We are now done with the vertical mass transports, so now h_new is h_sub
909  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
910  h_pre(i,j,k) = h_new(i,j,k)
911  enddo ; enddo ; enddo
912  call pass_var(h_pre,g%Domain)
913 
914  ! Second zonal and meridional advection
915  call update_h_horizontal_flux(g, gv, uhtr_sub, vhtr_sub, h_pre, h_new)
916  do k = 1, nz ; do i = is-1, ie+1 ; do j=js-1, je+1
917  h_vol(i,j,k) = h_pre(i,j,k)*g%areaT(i,j)
918  enddo ; enddo ; enddo
919  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, dt_iter, g, gv, &
920  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=30, x_first_in=x_before_y)
921 
922  ! Done with horizontal so now h_pre should be h_new
923  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
924  h_pre(i,j,k) = h_new(i,j,k)
925  enddo ; enddo ; enddo
926 
927  endif
928 
929  if (.not. z_first) then
930 
931  ! First zonal and meridional advection
932  call update_h_horizontal_flux(g, gv, uhtr_sub, vhtr_sub, h_pre, h_new)
933  do k = 1, nz ; do i = is-1, ie+1 ; do j=js-1, je+1
934  h_vol(i,j,k) = h_pre(i,j,k)*g%areaT(i,j)
935  enddo ; enddo ; enddo
936  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, dt_iter, g, gv, &
937  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=30, x_first_in=x_before_y)
938 
939  ! Done with horizontal so now h_pre should be h_new
940  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
941  h_pre(i,j,k) = h_new(i,j,k)
942  enddo ; enddo ; enddo
943 
944  ! Second vertical advection
945  call update_h_vertical_flux(g, gv, eatr_sub, ebtr_sub, h_pre, h_new)
946  call call_tracer_column_fns(h_pre, h_new, eatr_sub, ebtr_sub, &
947  fluxes, cs%mld, dt_iter, g, gv, cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
948  ! We are now done with the vertical mass transports, so now h_new is h_sub
949  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
950  h_pre(i,j,k) = h_new(i,j,k)
951  enddo ; enddo ; enddo
952 
953  endif
954 
955  ! Update remaining transports
956  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
957  eatr(i,j,k) = eatr(i,j,k) - eatr_sub(i,j,k)
958  ebtr(i,j,k) = ebtr(i,j,k) - ebtr_sub(i,j,k)
959  enddo ; enddo ; enddo
960 
961  do k = 1, nz ; do j=js-1,je+1 ; do i=is-2,ie+1
962  uhtr(i,j,k) = uhtr(i,j,k) - uhtr_sub(i,j,k)
963  enddo ; enddo ; enddo
964 
965  do k = 1, nz ; do j=js-2,je+1 ; do i=is-1,ie+1
966  vhtr(i,j,k) = vhtr(i,j,k) - vhtr_sub(i,j,k)
967  enddo ; enddo ; enddo
968 
969  call pass_var(eatr,g%Domain)
970  call pass_var(ebtr,g%Domain)
971  call pass_var(h_pre,g%Domain)
972  call pass_vector(uhtr,vhtr,g%Domain)
973  !
974  ! Calculate how close we are to converging by summing the remaining fluxes at each point
975  sum_abs_fluxes = 0.0
976  sum_u = 0.0
977  sum_v = 0.0
978  do k=1,nz; do j=js,je; do i=is,ie
979  sum_u = sum_u + abs(uhtr(i-1,j,k))+abs(uhtr(i,j,k))
980  sum_v = sum_v + abs(vhtr(i,j-1,k))+abs(vhtr(i,j,k))
981  sum_abs_fluxes = sum_abs_fluxes + abs(eatr(i,j,k)) + abs(ebtr(i,j,k)) + abs(uhtr(i-1,j,k)) + &
982  abs(uhtr(i,j,k)) + abs(vhtr(i,j-1,k)) + abs(vhtr(i,j,k))
983  enddo ; enddo ; enddo
984  call sum_across_pes(sum_abs_fluxes)
985 
986  write(mesg,*) "offline_advection_layer: Remaining u-flux, v-flux:", sum_u, sum_v
987  call mom_mesg(mesg)
988  if (sum_abs_fluxes==0) then
989  write(mesg,*) 'offline_advection_layer: Converged after iteration', iter
990  call mom_mesg(mesg)
991  exit
992  endif
993 
994  ! Switch order of Strang split every iteration
995  z_first = .not. z_first
996  x_before_y = .not. x_before_y
997  enddo
998 
999 end subroutine offline_advection_layer
1000 
1001 !> Update fields used in this round of offline transport. First fields are updated from files or from arrays
1002 !! read during initialization. Then if in an ALE-dependent coordinate, regrid/remap fields.
1003 subroutine update_offline_fields(CS, h, fluxes, do_ale)
1004  type(offline_transport_cs), pointer :: cs !< Control structure for offline module
1005  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: h !< The regridded layer thicknesses
1006  type(forcing), intent(inout) :: fluxes !< Pointers to forcing fields
1007  logical, intent(in ) :: do_ale !< True if using ALE
1008  ! Local variables
1009  integer :: i, j, k, is, ie, js, je, nz
1010  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: h_start
1011  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec ; nz = cs%GV%ke
1012 
1013  call cpu_clock_begin(cs%id_clock_read_fields)
1014  call calltree_enter("update_offline_fields, MOM_offline_main.F90")
1015 
1016  ! Store a copy of the layer thicknesses before ALE regrid/remap
1017  h_start(:,:,:) = h(:,:,:)
1018 
1019  ! Most fields will be read in from files
1020  call update_offline_from_files( cs%G, cs%GV, cs%nk_input, cs%mean_file, cs%sum_file, cs%snap_file, cs%surf_file, &
1021  cs%h_end, cs%uhtr, cs%vhtr, cs%tv%T, cs%tv%S, cs%mld, cs%Kd, fluxes, &
1022  cs%ridx_sum, cs%ridx_snap, cs%read_mld, cs%read_sw, .not. cs%read_all_ts_uvh, do_ale)
1023  ! If uh, vh, h_end, temp, salt were read in at the beginning, fields are copied from those arrays
1024  if (cs%read_all_ts_uvh) then
1025  call update_offline_from_arrays(cs%G, cs%GV, cs%nk_input, cs%ridx_sum, cs%mean_file, cs%sum_file, cs%snap_file, &
1026  cs%uhtr, cs%vhtr, cs%h_end, cs%uhtr_all, cs%vhtr_all, cs%hend_all, cs%tv%T, cs%tv%S, cs%temp_all, cs%salt_all)
1027  endif
1028  if (cs%debug) then
1029  call uvchksum("[uv]h after update offline from files and arrays", cs%uhtr, cs%vhtr, cs%G%HI)
1030  endif
1031 
1032  ! If using an ALE-dependent vertical coordinate, fields will need to be remapped
1033  if (do_ale) then
1034  ! These halo passes are necessary because u, v fields will need information 1 step into the halo
1035  call pass_var(h, cs%G%Domain)
1036  call pass_var(cs%tv%T, cs%G%Domain)
1037  call pass_var(cs%tv%S, cs%G%Domain)
1038  call ale_offline_inputs(cs%ALE_CSp, cs%G, cs%GV, h, cs%tv, cs%tracer_Reg, cs%uhtr, cs%vhtr, cs%Kd, cs%debug)
1039  if (cs%id_temp_regrid>0) call post_data(cs%id_temp_regrid, cs%tv%T, cs%diag)
1040  if (cs%id_salt_regrid>0) call post_data(cs%id_salt_regrid, cs%tv%S, cs%diag)
1041  if (cs%id_uhtr_regrid>0) call post_data(cs%id_uhtr_regrid, cs%uhtr, cs%diag)
1042  if (cs%id_vhtr_regrid>0) call post_data(cs%id_vhtr_regrid, cs%vhtr, cs%diag)
1043  if (cs%id_h_regrid>0) call post_data(cs%id_h_regrid, h, cs%diag)
1044  if (cs%debug) then
1045  call uvchksum("[uv]h after ALE regridding/remapping of inputs", cs%uhtr, cs%vhtr, cs%G%HI)
1046  call hchksum(h_start,"h_start after update offline from files and arrays", cs%G%HI)
1047  endif
1048  endif
1049 
1050  ! Update halos for some
1051  call pass_var(cs%h_end, cs%G%Domain)
1052  call pass_var(cs%tv%T, cs%G%Domain)
1053  call pass_var(cs%tv%S, cs%G%Domain)
1054 
1055  ! Update the read indices
1056  cs%ridx_snap = next_modulo_time(cs%ridx_snap,cs%numtime)
1057  cs%ridx_sum = next_modulo_time(cs%ridx_sum,cs%numtime)
1058 
1059  ! Apply masks/factors at T, U, and V points
1060  do k=1,nz ; do j=js,je ; do i=is,ie
1061  if (cs%G%mask2dT(i,j)<1.0) then
1062  cs%h_end(i,j,k) = cs%GV%Angstrom_H
1063  endif
1064  enddo ; enddo ; enddo
1065 
1066  do k=1,nz+1 ; do j=js,je ; do i=is,ie
1067  cs%Kd(i,j,k) = max(0.0, cs%Kd(i,j,k))
1068  if (cs%Kd_max>0.) then
1069  cs%Kd(i,j,k) = min(cs%Kd_max, cs%Kd(i,j,k))
1070  endif
1071  enddo ; enddo ; enddo
1072 
1073  do k=1,nz ; do j=js-1,je ; do i=is,ie
1074  if (cs%G%mask2dCv(i,j)<1.0) then
1075  cs%vhtr(i,j,k) = 0.0
1076  endif
1077  enddo ; enddo ; enddo
1078 
1079  do k=1,nz ; do j=js,je ; do i=is-1,ie
1080  if (cs%G%mask2dCu(i,j)<1.0) then
1081  cs%uhtr(i,j,k) = 0.0
1082  endif
1083  enddo ; enddo ; enddo
1084 
1085  if (cs%debug) then
1086  call uvchksum("[uv]htr_sub after update_offline_fields", cs%uhtr, cs%vhtr, cs%G%HI)
1087  call hchksum(cs%h_end, "h_end after update_offline_fields", cs%G%HI)
1088  call hchksum(cs%tv%T, "Temp after update_offline_fields", cs%G%HI)
1089  call hchksum(cs%tv%S, "Salt after update_offline_fields", cs%G%HI)
1090  endif
1091 
1092  call calltree_leave("update_offline_fields")
1093  call cpu_clock_end(cs%id_clock_read_fields)
1094 
1095 end subroutine update_offline_fields
1096 
1097 !> Initialize additional diagnostics required for offline tracer transport
1098 subroutine register_diags_offline_transport(Time, diag, CS)
1100  type(offline_transport_cs), pointer :: cs !< Control structure for offline module
1101  type(time_type), intent(in) :: time !< current model time
1102  type(diag_ctrl), intent(in) :: diag !< Structure that regulates diagnostic output
1103 
1104  ! U-cell fields
1105  cs%id_uhr = register_diag_field('ocean_model', 'uhr', diag%axesCuL, time, &
1106  'Zonal thickness fluxes remaining at end of advection', 'kg')
1107  cs%id_uhr_redist = register_diag_field('ocean_model', 'uhr_redist', diag%axesCuL, time, &
1108  'Zonal thickness fluxes to be redistributed vertically', 'kg')
1109  cs%id_uhr_end = register_diag_field('ocean_model', 'uhr_end', diag%axesCuL, time, &
1110  'Zonal thickness fluxes at end of offline step', 'kg')
1111 
1112  ! V-cell fields
1113  cs%id_vhr = register_diag_field('ocean_model', 'vhr', diag%axesCvL, time, &
1114  'Meridional thickness fluxes remaining at end of advection', 'kg')
1115  cs%id_vhr_redist = register_diag_field('ocean_model', 'vhr_redist', diag%axesCvL, time, &
1116  'Meridional thickness to be redistributed vertically', 'kg')
1117  cs%id_vhr_end = register_diag_field('ocean_model', 'vhr_end', diag%axesCvL, time, &
1118  'Meridional thickness at end of offline step', 'kg')
1119 
1120  ! T-cell fields
1121  cs%id_hdiff = register_diag_field('ocean_model', 'hdiff', diag%axesTL, time, &
1122  'Difference between the stored and calculated layer thickness', 'm')
1123  cs%id_hr = register_diag_field('ocean_model', 'hr', diag%axesTL, time, &
1124  'Layer thickness at end of offline step', 'm')
1125  cs%id_ear = register_diag_field('ocean_model', 'ear', diag%axesTL, time, &
1126  'Remaining thickness entrained from above', 'm')
1127  cs%id_ebr = register_diag_field('ocean_model', 'ebr', diag%axesTL, time, &
1128  'Remaining thickness entrained from below', 'm')
1129  cs%id_eta_pre_distribute = register_diag_field('ocean_model','eta_pre_distribute', &
1130  diag%axesT1, time, 'Total water column height before residual transport redistribution','m')
1131  cs%id_eta_post_distribute = register_diag_field('ocean_model','eta_post_distribute', &
1132  diag%axesT1, time, 'Total water column height after residual transport redistribution','m')
1133  cs%id_eta_diff_end = register_diag_field('ocean_model','eta_diff_end', diag%axesT1, time, &
1134  'Difference in total water column height from online and offline ' // &
1135  'at the end of the offline timestep','m')
1136  cs%id_h_redist = register_diag_field('ocean_model','h_redist', diag%axesTL, time, &
1137  'Layer thicknesses before redistribution of mass fluxes','m')
1138 
1139  ! Regridded/remapped input fields
1140  cs%id_uhtr_regrid = register_diag_field('ocean_model', 'uhtr_regrid', diag%axesCuL, time, &
1141  'Zonal mass transport regridded/remapped onto offline grid','kg')
1142  cs%id_vhtr_regrid = register_diag_field('ocean_model', 'vhtr_regrid', diag%axesCvL, time, &
1143  'Meridional mass transport regridded/remapped onto offline grid','kg')
1144  cs%id_temp_regrid = register_diag_field('ocean_model', 'temp_regrid', diag%axesTL, time, &
1145  'Temperature regridded/remapped onto offline grid','C')
1146  cs%id_salt_regrid = register_diag_field('ocean_model', 'salt_regrid', diag%axesTL, time, &
1147  'Salinity regridded/remapped onto offline grid','g kg-1')
1148  cs%id_h_regrid = register_diag_field('ocean_model', 'h_regrid', diag%axesTL, time, &
1149  'Layer thicknesses regridded/remapped onto offline grid','m')
1150 
1151 
1152 end subroutine register_diags_offline_transport
1153 
1154 !> Posts diagnostics related to offline convergence diagnostics
1155 subroutine post_offline_convergence_diags(CS, h_off, h_end, uhtr, vhtr)
1156  type(offline_transport_cs), intent(in ) :: cs !< Offline control structure
1157  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_off !< Thicknesses at end of offline step
1158  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_end !< Stored thicknesses
1159  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: uhtr !< Remaining zonal mass transport
1160  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(inout) :: vhtr !< Remaining meridional mass transport
1161 
1162  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_diff
1163  integer :: i, j, k
1164 
1165  if (cs%id_eta_diff_end>0) then
1166  ! Calculate difference in column thickness
1167  eta_diff = 0.
1168  do k=1,cs%GV%ke ; do j=cs%G%jsc,cs%G%jec ; do i=cs%G%isc,cs%G%iec
1169  eta_diff(i,j) = eta_diff(i,j) + h_off(i,j,k)
1170  enddo ; enddo ; enddo
1171  do k=1,cs%GV%ke ; do j=cs%G%jsc,cs%G%jec ; do i=cs%G%isc,cs%G%iec
1172  eta_diff(i,j) = eta_diff(i,j) - h_end(i,j,k)
1173  enddo ; enddo ; enddo
1174 
1175  call post_data(cs%id_eta_diff_end, eta_diff, cs%diag)
1176  endif
1177 
1178  if (cs%id_hdiff>0) call post_data(cs%id_hdiff, h_off-h_end, cs%diag)
1179  if (cs%id_hr>0) call post_data(cs%id_hr, h_off, cs%diag)
1180  if (cs%id_uhr_end>0) call post_data(cs%id_uhr_end, uhtr, cs%diag)
1181  if (cs%id_vhr_end>0) call post_data(cs%id_vhr_end, vhtr, cs%diag)
1182 
1183 end subroutine post_offline_convergence_diags
1184 
1185 !> Extracts members of the offline main control structure. All arguments are optional except
1186 !! the control structure itself
1187 subroutine extract_offline_main(CS, uhtr, vhtr, eatr, ebtr, h_end, accumulated_time, &
1188  dt_offline, dt_offline_vertical, skip_diffusion)
1189  type(offline_transport_cs), target, intent(in ) :: cs !< Offline control structure
1190  ! Returned optional arguments
1191  real, dimension(:,:,:), optional, pointer :: uhtr !< Remaining zonal mass transport [H m2 ~> m3 or kg]
1192  real, dimension(:,:,:), optional, pointer :: vhtr !< Remaining meridional mass transport [H m2 ~> m3 or kg]
1193  real, dimension(:,:,:), optional, pointer :: eatr !< Amount of fluid entrained from the layer above within
1194  !! one time step [H ~> m or kg m-2]
1195  real, dimension(:,:,:), optional, pointer :: ebtr !< Amount of fluid entrained from the layer below within
1196  !! one time step [H ~> m or kg m-2]
1197  real, dimension(:,:,:), optional, pointer :: h_end !< Thicknesses at the end of offline timestep
1198  !! [H ~> m or kg m-2]
1199  !### Why are the following variables integers?
1200  integer, optional, pointer :: accumulated_time !< Length of time accumulated in the
1201  !! current offline interval [s]
1202  integer, optional, intent( out) :: dt_offline !< Timestep used for offline tracers [s]
1203  integer, optional, intent( out) :: dt_offline_vertical !< Timestep used for calls to tracer
1204  !! vertical physics [s]
1205  logical, optional, intent( out) :: skip_diffusion !< Skips horizontal diffusion of tracers
1206 
1207  ! Pointers to 3d members
1208  if (present(uhtr)) uhtr => cs%uhtr
1209  if (present(vhtr)) vhtr => cs%vhtr
1210  if (present(eatr)) eatr => cs%eatr
1211  if (present(ebtr)) ebtr => cs%ebtr
1212  if (present(h_end)) h_end => cs%h_end
1213 
1214  ! Pointers to integer members which need to be modified
1215  if (present(accumulated_time)) accumulated_time => cs%accumulated_time
1216 
1217  ! Return value of non-modified integers
1218  if (present(dt_offline)) dt_offline = cs%dt_offline
1219  if (present(dt_offline_vertical)) dt_offline_vertical = cs%dt_offline_vertical
1220  if (present(skip_diffusion)) skip_diffusion = cs%skip_diffusion
1221 
1222 end subroutine extract_offline_main
1223 
1224 !> Inserts (assigns values to) members of the offline main control structure. All arguments
1225 !! are optional except for the CS itself
1226 subroutine insert_offline_main(CS, ALE_CSp, diabatic_CSp, diag, OBC, tracer_adv_CSp, &
1227  tracer_flow_CSp, tracer_Reg, tv, G, GV, x_before_y, debug)
1228  type(offline_transport_cs), intent(inout) :: cs !< Offline control structure
1229  ! Inserted optional arguments
1230  type(ale_cs), &
1231  target, optional, intent(in ) :: ale_csp !< A pointer to the ALE control structure
1232  type(diabatic_cs), &
1233  target, optional, intent(in ) :: diabatic_csp !< A pointer to the diabatic control structure
1234  type(diag_ctrl), &
1235  target, optional, intent(in ) :: diag !< A pointer to the structure that regulates diagnostic output
1236  type(ocean_obc_type), &
1237  target, optional, intent(in ) :: obc !< A pointer to the open boundary condition control structure
1238  type(tracer_advect_cs), &
1239  target, optional, intent(in ) :: tracer_adv_csp !< A pointer to the tracer advection control structure
1240  type(tracer_flow_control_cs), &
1241  target, optional, intent(in ) :: tracer_flow_csp !< A pointer to the tracer flow control control structure
1242  type(tracer_registry_type), &
1243  target, optional, intent(in ) :: tracer_reg !< A pointer to the tracer registry
1244  type(thermo_var_ptrs), &
1245  target, optional, intent(in ) :: tv !< A structure pointing to various thermodynamic variables
1246  type(ocean_grid_type), &
1247  target, optional, intent(in ) :: g !< ocean grid structure
1248  type(verticalgrid_type), &
1249  target, optional, intent(in ) :: gv !< ocean vertical grid structure
1250  logical, optional, intent(in ) :: x_before_y !< Indicates which horizontal direction is advected first
1251  logical, optional, intent(in ) :: debug !< If true, write verbose debugging messages
1252 
1253 
1254  if (present(ale_csp)) cs%ALE_CSp => ale_csp
1255  if (present(diabatic_csp)) cs%diabatic_CSp => diabatic_csp
1256  if (present(diag)) cs%diag => diag
1257  if (present(obc)) cs%OBC => obc
1258  if (present(tracer_adv_csp)) cs%tracer_adv_CSp => tracer_adv_csp
1259  if (present(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
1260  if (present(tracer_reg)) cs%tracer_Reg => tracer_reg
1261  if (present(tv)) cs%tv => tv
1262  if (present(g)) cs%G => g
1263  if (present(gv)) cs%GV => gv
1264  if (present(x_before_y)) cs%x_before_y = x_before_y
1265  if (present(debug)) cs%debug = debug
1266 
1267 end subroutine insert_offline_main
1268 
1269 !> Initializes the control structure for offline transport and reads in some of the
1270 ! run time parameters from MOM_input
1271 subroutine offline_transport_init(param_file, CS, diabatic_CSp, G, GV)
1273  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1274  type(offline_transport_cs), pointer :: cs !< Offline control structure
1275  type(diabatic_cs), intent(in) :: diabatic_csp !< The diabatic control structure
1276  type(ocean_grid_type), target, intent(in) :: g !< ocean grid structure
1277  type(verticalgrid_type), target, intent(in) :: gv !< ocean vertical grid structure
1278 
1279  character(len=40) :: mdl = "offline_transport"
1280  character(len=20) :: redistribute_method
1281 
1282  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
1283  integer :: isdb, iedb, jsdb, jedb
1284 
1285  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
1286  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1287  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1288 
1289  call calltree_enter("offline_transport_init, MOM_offline_control.F90")
1290 
1291  if (associated(cs)) then
1292  call mom_error(warning, "offline_transport_init called with an associated "// &
1293  "control structure.")
1294  return
1295  endif
1296  allocate(cs)
1297  call log_version(param_file, mdl,version, "This module allows for tracers to be run offline")
1298 
1299  ! Parse MOM_input for offline control
1300  call get_param(param_file, mdl, "OFFLINEDIR", cs%offlinedir, &
1301  "Input directory where the offline fields can be found", fail_if_missing = .true.)
1302  call get_param(param_file, mdl, "OFF_SUM_FILE", cs%sum_file, &
1303  "Filename where the accumulated fields can be found", fail_if_missing = .true.)
1304  call get_param(param_file, mdl, "OFF_SNAP_FILE", cs%snap_file, &
1305  "Filename where snapshot fields can be found", fail_if_missing = .true.)
1306  call get_param(param_file, mdl, "OFF_MEAN_FILE", cs%mean_file, &
1307  "Filename where averaged fields can be found", fail_if_missing = .true.)
1308  call get_param(param_file, mdl, "OFF_SURF_FILE", cs%surf_file, &
1309  "Filename where averaged fields can be found", fail_if_missing = .true.)
1310  call get_param(param_file, mdl, "NUMTIME", cs%numtime, &
1311  "Number of timelevels in offline input files", fail_if_missing = .true.)
1312  call get_param(param_file, mdl, "NK_INPUT", cs%nk_input, &
1313  "Number of vertical levels in offline input files", default = nz)
1314  call get_param(param_file, mdl, "DT_OFFLINE", cs%dt_offline, &
1315  "Length of time between reading in of input fields", fail_if_missing = .true.)
1316  call get_param(param_file, mdl, "DT_OFFLINE_VERTICAL", cs%dt_offline_vertical, &
1317  "Length of the offline timestep for tracer column sources/sinks " //&
1318  "This should be set to the length of the coupling timestep for " //&
1319  "tracers which need shortwave fluxes", fail_if_missing = .true.)
1320  call get_param(param_file, mdl, "START_INDEX", cs%start_index, &
1321  "Which time index to start from", default=1)
1322  call get_param(param_file, mdl, "FIELDS_ARE_OFFSET", cs%fields_are_offset, &
1323  "True if the time-averaged fields and snapshot fields "//&
1324  "are offset by one time level", default=.false.)
1325  call get_param(param_file, mdl, "REDISTRIBUTE_METHOD", redistribute_method, &
1326  "Redistributes any remaining horizontal fluxes throughout " //&
1327  "the rest of water column. Options are 'barotropic' which " //&
1328  "evenly distributes flux throughout the entire water column, " //&
1329  "'upwards' which adds the maximum of the remaining flux in " //&
1330  "each layer above, both which first applies upwards and then " //&
1331  "barotropic, and 'none' which does no redistribution", &
1332  default='barotropic')
1333  call get_param(param_file, mdl, "NUM_OFF_ITER", cs%num_off_iter, &
1334  "Number of iterations to subdivide the offline tracer advection and diffusion", &
1335  default = 60)
1336  call get_param(param_file, mdl, "OFF_ALE_MOD", cs%off_ale_mod, &
1337  "Sets how many horizontal advection steps are taken before an ALE " //&
1338  "remapping step is done. 1 would be x->y->ALE, 2 would be" //&
1339  "x->y->x->y->ALE", default = 1)
1340  call get_param(param_file, mdl, "PRINT_ADV_OFFLINE", cs%print_adv_offline, &
1341  "Print diagnostic output every advection subiteration",default=.false.)
1342  call get_param(param_file, mdl, "SKIP_DIFFUSION_OFFLINE", cs%skip_diffusion, &
1343  "Do not do horizontal diffusion",default=.false.)
1344  call get_param(param_file, mdl, "READ_SW", cs%read_sw, &
1345  "Read in shortwave radiation field instead of using values from the coupler"//&
1346  "when in offline tracer mode",default=.false.)
1347  call get_param(param_file, mdl, "READ_MLD", cs%read_mld, &
1348  "Read in mixed layer depths for tracers which exchange with the atmosphere"//&
1349  "when in offline tracer mode",default=.false.)
1350  call get_param(param_file, mdl, "MLD_VAR_NAME", cs%mld_var_name, &
1351  "Name of the variable containing the depth of active mixing",&
1352  default='ePBL_h_ML')
1353  call get_param(param_file, mdl, "OFFLINE_ADD_DIURNAL_SW", cs%diurnal_sw, &
1354  "Adds a synthetic diurnal cycle in the same way that the ice " // &
1355  "model would have when time-averaged fields of shortwave " // &
1356  "radiation are read in", default=.false.)
1357  call get_param(param_file, mdl, "KD_MAX", cs%Kd_max, &
1358  "The maximum permitted increment for the diapycnal "//&
1359  "diffusivity from TKE-based parameterizations, or a "//&
1360  "negative value for no limit.", units="m2 s-1", default=-1.0)
1361  call get_param(param_file, mdl, "MIN_RESIDUAL_TRANSPORT", cs%min_residual, &
1362  "How much remaining transport before the main offline advection "// &
1363  "is exited. The default value corresponds to about 1 meter of " // &
1364  "difference in a grid cell", default = 1.e9)
1365  call get_param(param_file, mdl, "READ_ALL_TS_UVH", cs%read_all_ts_uvh, &
1366  "Reads all time levels of a subset of the fields necessary to run " // &
1367  "the model offline. This can require a large amount of memory "// &
1368  "and will make initialization very slow. However, for offline "// &
1369  "runs spanning more than a year this can reduce total I/O overhead", &
1370  default = .false.)
1371 
1372  ! Concatenate offline directory and file names
1373  cs%snap_file = trim(cs%offlinedir)//trim(cs%snap_file)
1374  cs%mean_file = trim(cs%offlinedir)//trim(cs%mean_file)
1375  cs%sum_file = trim(cs%offlinedir)//trim(cs%sum_file)
1376  cs%surf_file = trim(cs%offlinedir)//trim(cs%surf_file)
1377 
1378  cs%num_vert_iter = cs%dt_offline/cs%dt_offline_vertical
1379 
1380  ! Map redistribute_method onto logicals in CS
1381  select case (redistribute_method)
1382  case ('barotropic')
1383  cs%redistribute_barotropic = .true.
1384  cs%redistribute_upwards = .false.
1385  case ('upwards')
1386  cs%redistribute_barotropic = .false.
1387  cs%redistribute_upwards = .true.
1388  case ('both')
1389  cs%redistribute_barotropic = .true.
1390  cs%redistribute_upwards = .true.
1391  case ('none')
1392  cs%redistribute_barotropic = .false.
1393  cs%redistribute_upwards = .false.
1394  end select
1395 
1396  ! Set the accumulated time to zero
1397  cs%accumulated_time = 0
1398  ! Set the starting read index for time-averaged and snapshotted fields
1399  cs%ridx_sum = cs%start_index
1400  if (cs%fields_are_offset) cs%ridx_snap = next_modulo_time(cs%start_index,cs%numtime)
1401  if (.not. cs%fields_are_offset) cs%ridx_snap = cs%start_index
1402 
1403  ! Copy members from other modules
1404  call extract_diabatic_member(diabatic_csp, opacity_csp=cs%opacity_CSp, optics_csp=cs%optics, &
1405  diabatic_aux_csp=cs%diabatic_aux_CSp, &
1406  evap_cfl_limit=cs%evap_CFL_limit, &
1407  minimum_forcing_depth=cs%minimum_forcing_depth)
1408 
1409  ! Grid pointer assignments
1410  cs%G => g
1411  cs%GV => gv
1412 
1413  ! Allocate arrays
1414  allocate(cs%uhtr(isdb:iedb,jsd:jed,nz)) ; cs%uhtr(:,:,:) = 0.0
1415  allocate(cs%vhtr(isd:ied,jsdb:jedb,nz)) ; cs%vhtr(:,:,:) = 0.0
1416  allocate(cs%eatr(isd:ied,jsd:jed,nz)) ; cs%eatr(:,:,:) = 0.0
1417  allocate(cs%ebtr(isd:ied,jsd:jed,nz)) ; cs%ebtr(:,:,:) = 0.0
1418  allocate(cs%h_end(isd:ied,jsd:jed,nz)) ; cs%h_end(:,:,:) = 0.0
1419  allocate(cs%netMassOut(g%isd:g%ied,g%jsd:g%jed)) ; cs%netMassOut(:,:) = 0.0
1420  allocate(cs%netMassIn(g%isd:g%ied,g%jsd:g%jed)) ; cs%netMassIn(:,:) = 0.0
1421  allocate(cs%Kd(isd:ied,jsd:jed,nz+1)) ; cs%Kd = 0.
1422  if (cs%read_mld) then
1423  allocate(cs%mld(g%isd:g%ied,g%jsd:g%jed)) ; cs%mld(:,:) = 0.0
1424  endif
1425 
1426  if (cs%read_all_ts_uvh) then
1427  call read_all_input(cs)
1428  endif
1429 
1430  ! Initialize ids for clocks used in offline routines
1431  cs%id_clock_read_fields = cpu_clock_id('(Offline read fields)',grain=clock_module)
1432  cs%id_clock_offline_diabatic = cpu_clock_id('(Offline diabatic)',grain=clock_module)
1433  cs%id_clock_offline_adv = cpu_clock_id('(Offline transport)',grain=clock_module)
1434  cs%id_clock_redistribute = cpu_clock_id('(Offline redistribute)',grain=clock_module)
1435 
1436  call calltree_leave("offline_transport_init")
1437 
1438 end subroutine offline_transport_init
1439 
1440 !> Coordinates the allocation and reading in all time levels of uh, vh, hend, temp, and salt from files. Used
1441 !! when read_all_ts_uvh
1442 subroutine read_all_input(CS)
1443  type(offline_transport_cs), intent(inout) :: CS !< Control structure for offline module
1444 
1445  integer :: is, ie, js, je, isd, ied, jsd, jed, nz, t, ntime
1446  integer :: IsdB, IedB, JsdB, JedB
1447 
1448  nz = cs%GV%ke ; ntime = cs%numtime
1449  isd = cs%G%isd ; ied = cs%G%ied ; jsd = cs%G%jsd ; jed = cs%G%jed
1450  isdb = cs%G%IsdB ; iedb = cs%G%IedB ; jsdb = cs%G%JsdB ; jedb = cs%G%JedB
1451 
1452  ! Extra safety check that we're not going to overallocate any arrays
1453  if (cs%read_all_ts_uvh) then
1454  if (allocated(cs%uhtr_all)) call mom_error(fatal, "uhtr_all is already allocated")
1455  if (allocated(cs%vhtr_all)) call mom_error(fatal, "vhtr_all is already allocated")
1456  if (allocated(cs%hend_all)) call mom_error(fatal, "hend_all is already allocated")
1457  if (allocated(cs%temp_all)) call mom_error(fatal, "temp_all is already allocated")
1458  if (allocated(cs%salt_all)) call mom_error(fatal, "salt_all is already allocated")
1459 
1460  allocate(cs%uhtr_all(isdb:iedb,jsd:jed,nz,ntime)) ; cs%uhtr_all(:,:,:,:) = 0.0
1461  allocate(cs%vhtr_all(isd:ied,jsdb:jedb,nz,ntime)) ; cs%vhtr_all(:,:,:,:) = 0.0
1462  allocate(cs%hend_all(isd:ied,jsd:jed,nz,ntime)) ; cs%hend_all(:,:,:,:) = 0.0
1463  allocate(cs%temp_all(isd:ied,jsd:jed,nz,1:ntime)) ; cs%temp_all(:,:,:,:) = 0.0
1464  allocate(cs%salt_all(isd:ied,jsd:jed,nz,1:ntime)) ; cs%salt_all(:,:,:,:) = 0.0
1465 
1466  call mom_mesg("Reading in uhtr, vhtr, h_start, h_end, temp, salt")
1467  do t = 1,ntime
1468  call mom_read_vector(cs%snap_file, 'uhtr_sum', 'vhtr_sum', cs%uhtr_all(:,:,1:cs%nk_input,t), &
1469  cs%vhtr_all(:,:,1:cs%nk_input,t), cs%G%Domain, timelevel=t)
1470  call mom_read_data(cs%snap_file,'h_end', cs%hend_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1471  timelevel=t, position=center)
1472  call mom_read_data(cs%mean_file,'temp', cs%temp_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1473  timelevel=t, position=center)
1474  call mom_read_data(cs%mean_file,'salt', cs%salt_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1475  timelevel=t, position=center)
1476  enddo
1477  endif
1478 
1479 end subroutine read_all_input
1480 
1481 !> Deallocates (if necessary) arrays within the offline control structure
1482 subroutine offline_transport_end(CS)
1483  type(offline_transport_cs), pointer :: cs !< Control structure for offline module
1484 
1485  ! Explicitly allocate all allocatable arrays
1486  deallocate(cs%uhtr)
1487  deallocate(cs%vhtr)
1488  deallocate(cs%eatr)
1489  deallocate(cs%ebtr)
1490  deallocate(cs%h_end)
1491  deallocate(cs%netMassOut)
1492  deallocate(cs%netMassIn)
1493  deallocate(cs%Kd)
1494  if (cs%read_mld) deallocate(cs%mld)
1495  if (cs%read_all_ts_uvh) then
1496  deallocate(cs%uhtr_all)
1497  deallocate(cs%vhtr_all)
1498  deallocate(cs%hend_all)
1499  deallocate(cs%temp_all)
1500  deallocate(cs%salt_all)
1501  endif
1502 
1503  deallocate(cs)
1504 
1505 end subroutine offline_transport_end
1506 
1507 !> \namespace mom_offline_main
1508 !! \section offline_overview Offline Tracer Transport in MOM6
1509 !! 'Offline tracer modeling' uses physical fields (e.g. mass transports and layer thicknesses) saved
1510 !! from a previous integration of the physical model to transport passive tracers. These fields are
1511 !! accumulated or averaged over a period of time (in this test case, 1 day) and used to integrate
1512 !! portions of the MOM6 code base that handle the 3d advection and diffusion of passive tracers.
1513 !!
1514 !! The distribution of tracers in the ocean modeled offline should not be expected to match an online
1515 !! simulation. Accumulating transports over more than one online model timestep implicitly assumes
1516 !! homogeneity over that time period and essentially aliases over processes that occur with higher
1517 !! frequency. For example, consider the case of a surface boundary layer with a strong diurnal cycle.
1518 !! An offline simulation with a 1 day timestep, captures the net transport into or out of that layer,
1519 !! but not the exact cycling. This effective aliasing may also complicate online model configurations
1520 !! which strongly-eddying regions. In this case, the offline model timestep must be limited to some
1521 !! fraction of the eddy correlation timescale. Lastly, the nonlinear advection scheme which applies
1522 !! limited mass-transports over a sequence of iterations means that tracers are not transported along
1523 !! exactly the same path as they are in the online model.
1524 !!
1525 !! This capability has currently targeted the Baltic_ALE_z test case, though some work has also been
1526 !! done with the OM4 1/2 degree configuration. Work is ongoing to develop recommendations and best
1527 !! practices for investigators seeking to use MOM6 for offline tracer modeling.
1528 !!
1529 !! \section offline_technical Implementation of offline routine in MOM6
1530 !!
1531 !! The subroutine step_tracers that coordinates this can be found in MOM.F90 and is only called
1532 !! using the solo ocean driver. This is to avoid issues with coupling to other climate components
1533 !! that may be relying on fluxes from the ocean to be coupled more often than the offline time step.
1534 !! Other routines related to offline tracer modeling can be found in tracers/MOM_offline_control.F90
1535 !!
1536 !! As can also be seen in the comments for the step_tracers subroutine, an offline time step
1537 !! comprises the following steps:
1538 !! -# Using the layer thicknesses and tracer concentrations from the previous timestep,
1539 !! half of the accumulated vertical mixing (eatr and ebtr) is applied in the call to
1540 !! tracer_column_fns.
1541 !! For tracers whose source/sink terms need dt, this value is set to 1/2 dt_offline
1542 !! -# Half of the accumulated surface freshwater fluxes are applied
1543 !! START ITERATION
1544 !! -# Accumulated mass fluxes are used to do horizontal transport. The number of iterations
1545 !! used in advect_tracer is limited to 2 (e.g x->y->x->y). The remaining mass fluxes are
1546 !! stored for later use and resulting layer thicknesses fed into the next step
1547 !! -# Tracers and the h-grid are regridded and remapped in a call to ALE. This allows for
1548 !! layers which might 'vanish' because of horizontal mass transport to be 'reinflated'
1549 !! and essentially allows for the vertical transport of tracers
1550 !! -# Check that transport is done if the remaining mass fluxes equals 0 or if the max
1551 !! number of iterations has been reached
1552 !! END ITERATION
1553 !! -# Repeat steps 1 and 2
1554 !! -# Redistribute any residual mass fluxes that remain after the advection iterations
1555 !! in a barotropic manner, progressively upward through the water column.
1556 !! -# Force a remapping to the stored layer thicknesses that correspond to the snapshot of
1557 !! the online model at the end of an accumulation interval
1558 !! -# Reset T/S and h to their stored snapshotted values to prevent model drift
1559 !!
1560 !! \section offline_evaluation Evaluating the utility of an offline tracer model
1561 !! How well an offline tracer model can be used as an alternative to integrating tracers online
1562 !! with the prognostic model must be evaluated for each application. This efficacy may be related
1563 !! to the native coordinate of the online model, to the length of the offline timestep, and to the
1564 !! behavior of the tracer itself.
1565 !!
1566 !! A framework for formally regression testing the offline capability still needs to be developed.
1567 !! However, as a simple way of testing whether the offline model is nominally behaving as expected,
1568 !! the total inventory of the advection test tracers (tr1, tr2, etc.) should be conserved between
1569 !! time steps except for the last 4 decimal places. As a general guideline, an offline timestep of
1570 !! 5 days or less.
1571 !!
1572 !! \section offline_parameters Runtime parameters for offline tracers
1573 !! - OFFLINEDIR: Input directory where the offline fields can be found
1574 !! - OFF_SUM_FILE: Filename where the accumulated fields can be found (e.g. horizontal mass transports)
1575 !! - OFF_SNAP_FILE: Filename where snapshot fields can be found (e.g. end of timestep layer thickness)
1576 !! - START_INDEX: Which timelevel of the input files to read first
1577 !! - NUMTIME: How many timelevels to read before 'looping' back to 1
1578 !! - FIELDS_ARE_OFFSET: True if the time-averaged fields and snapshot fields are offset by one
1579 !! time level, probably not needed
1580 !! -NUM_OFF_ITER: Maximum number of iterations to do for the nonlinear advection scheme
1581 !! -REDISTRIBUTE_METHOD: Redistributes any remaining horizontal fluxes throughout the rest of water column.
1582 !! Options are 'barotropic' which "evenly distributes flux throughout the entire water
1583 !! column,'upwards' which adds the maximum of the remaining flux in each layer above,
1584 !! and 'none' which does no redistribution"
1585 
1586 end module mom_offline_main
1587 
mom_tracer_advect::tracer_advect_cs
Control structure for this module.
Definition: MOM_tracer_advect.F90:30
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_opacity::opacity_cs
The control structure with paramters for the MOM_opacity module.
Definition: MOM_opacity.F90:50
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_checksums::uvchksum
Checksums a pair velocity arrays (2d or 3d) staggered at C-grid locations.
Definition: MOM_checksums.F90:28
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_tracer_diabatic
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
Definition: MOM_tracer_diabatic.F90:4
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_checksums
Routines to calculate checksums of various array and vector types.
Definition: MOM_checksums.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_ale
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_io::mom_read_vector
Read a pair of data fields representing the two components of a vector from a file.
Definition: MOM_io.F90:82
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
mom_cpu_clock
Wraps the MPP cpu clock functions.
Definition: MOM_cpu_clock.F90:2
mom_offline_main
The routines here implement the offline tracer algorithm used in MOM6. These are called from step_off...
Definition: MOM_offline_main.F90:3
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_ale::ale_cs
ALE control structure.
Definition: MOM_ALE.F90:63
mom_tracer_flow_control
Orchestrates the registration and calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:2
mom_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
mom_opacity::optics_type
This type is used to store information about ocean optical properties.
Definition: MOM_opacity.F90:25
mom_diabatic_driver
This routine drives the diabatic/dianeutral physics for MOM.
Definition: MOM_diabatic_driver.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:186
mom_tracer_flow_control::tracer_flow_control_cs
The control structure for orchestrating the calling of tracer packages.
Definition: MOM_tracer_flow_control.F90:75
mom_offline_aux
Contains routines related to offline transport of tracers. These routines are likely to be called fro...
Definition: MOM_offline_aux.F90:3
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:49
mom_diabatic_aux::diabatic_aux_cs
Control structure for diabatic_aux.
Definition: MOM_diabatic_aux.F90:42
mom_tracer_advect
This program contains the subroutines that advect tracers along coordinate surfaces.
Definition: MOM_tracer_advect.F90:3
mom_offline_main::offline_transport_cs
The control structure for the offline transport module.
Definition: MOM_offline_main.F90:45
mom_checksums::hchksum
Checksums an array (2d or 3d) staggered at tracer points.
Definition: MOM_checksums.F90:48
mom_opacity
Routines used to calculate the opacity of the ocean.
Definition: MOM_opacity.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_diabatic_aux
Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surfa...
Definition: MOM_diabatic_aux.F90:3
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diabatic_driver::diabatic_cs
Control structure for this module.
Definition: MOM_diabatic_driver.F90:92
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239
mom_file_parser::read_param
An overloaded interface to read various types of parameters.
Definition: MOM_file_parser.F90:90