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