MOM6
mom_offline_main Module Reference

Detailed Description

The routines here implement the offline tracer algorithm used in MOM6. These are called from step_offline Some routines called here can be found in the MOM_offline_aux module.

Offline Tracer Transport in MOM6

'Offline tracer modeling' uses physical fields (e.g. mass transports and layer thicknesses) saved from a previous integration of the physical model to transport passive tracers. These fields are accumulated or averaged over a period of time (in this test case, 1 day) and used to integrate portions of the MOM6 code base that handle the 3d advection and diffusion of passive tracers.

The distribution of tracers in the ocean modeled offline should not be expected to match an online simulation. Accumulating transports over more than one online model timestep implicitly assumes homogeneity over that time period and essentially aliases over processes that occur with higher frequency. For example, consider the case of a surface boundary layer with a strong diurnal cycle. An offline simulation with a 1 day timestep, captures the net transport into or out of that layer, but not the exact cycling. This effective aliasing may also complicate online model configurations which strongly-eddying regions. In this case, the offline model timestep must be limited to some fraction of the eddy correlation timescale. Lastly, the nonlinear advection scheme which applies limited mass-transports over a sequence of iterations means that tracers are not transported along exactly the same path as they are in the online model.

This capability has currently targeted the Baltic_ALE_z test case, though some work has also been done with the OM4 1/2 degree configuration. Work is ongoing to develop recommendations and best practices for investigators seeking to use MOM6 for offline tracer modeling.

Implementation of offline routine in MOM6

The subroutine step_tracers that coordinates this can be found in MOM.F90 and is only called using the solo ocean driver. This is to avoid issues with coupling to other climate components that may be relying on fluxes from the ocean to be coupled more often than the offline time step. Other routines related to offline tracer modeling can be found in tracers/MOM_offline_control.F90

As can also be seen in the comments for the step_tracers subroutine, an offline time step comprises the following steps:

  1. Using the layer thicknesses and tracer concentrations from the previous timestep, half of the accumulated vertical mixing (eatr and ebtr) is applied in the call to tracer_column_fns. For tracers whose source/sink terms need dt, this value is set to 1/2 dt_offline
  2. Half of the accumulated surface freshwater fluxes are applied START ITERATION
  3. Accumulated mass fluxes are used to do horizontal transport. The number of iterations used in advect_tracer is limited to 2 (e.g x->y->x->y). The remaining mass fluxes are stored for later use and resulting layer thicknesses fed into the next step
  4. Tracers and the h-grid are regridded and remapped in a call to ALE. This allows for layers which might 'vanish' because of horizontal mass transport to be 'reinflated' and essentially allows for the vertical transport of tracers
  5. Check that transport is done if the remaining mass fluxes equals 0 or if the max number of iterations has been reached END ITERATION
  6. Repeat steps 1 and 2
  7. Redistribute any residual mass fluxes that remain after the advection iterations in a barotropic manner, progressively upward through the water column.
  8. Force a remapping to the stored layer thicknesses that correspond to the snapshot of the online model at the end of an accumulation interval
  9. Reset T/S and h to their stored snapshotted values to prevent model drift

Evaluating the utility of an offline tracer model

How well an offline tracer model can be used as an alternative to integrating tracers online with the prognostic model must be evaluated for each application. This efficacy may be related to the native coordinate of the online model, to the length of the offline timestep, and to the behavior of the tracer itself.

A framework for formally regression testing the offline capability still needs to be developed. However, as a simple way of testing whether the offline model is nominally behaving as expected, the total inventory of the advection test tracers (tr1, tr2, etc.) should be conserved between time steps except for the last 4 decimal places. As a general guideline, an offline timestep of 5 days or less.

Runtime parameters for offline tracers

  • OFFLINEDIR: Input directory where the offline fields can be found
  • OFF_SUM_FILE: Filename where the accumulated fields can be found (e.g. horizontal mass transports)
  • OFF_SNAP_FILE: Filename where snapshot fields can be found (e.g. end of timestep layer thickness)
  • START_INDEX: Which timelevel of the input files to read first
  • NUMTIME: How many timelevels to read before 'looping' back to 1
  • FIELDS_ARE_OFFSET: True if the time-averaged fields and snapshot fields are offset by one time level, probably not needed -NUM_OFF_ITER: Maximum number of iterations to do for the nonlinear advection scheme -REDISTRIBUTE_METHOD: Redistributes any remaining horizontal fluxes throughout the rest of water column. Options are 'barotropic' which "evenly distributes flux throughout the entire water column,'upwards' which adds the maximum of the remaining flux in each layer above, and 'none' which does no redistribution"

Data Types

type  offline_transport_cs
 The control structure for the offline transport module. More...
 

Functions/Subroutines

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 regridding/remapping step. The loop in this routine is exited if remaining residual transports are below a runtime-specified value or a maximum number of iterations is reached. More...
 
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 remaining transport. Two different ways are offered, 'barotropic' means that the residual is distributed equally throughout the water column. 'upwards' attempts to redistribute the transport in the layers above and will eventually work down the entire water column. More...
 
real function remaining_transport_sum (CS, uhtr, vhtr)
 Sums any non-negligible remaining transport to check for advection convergence. More...
 
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 interpolated vertical diffusivities are calculated and then any tracer column functions are done which can include vertical diffuvities and source/sink terms. More...
 
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 of the ocean) fluxes. More...
 
subroutine, public offline_fw_fluxes_out_ocean (G, GV, CS, fluxes, h, out_flux_optional)
 Apply negative freshwater fluxes (out of the ocean) More...
 
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. More...
 
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 arrays read during initialization. Then if in an ALE-dependent coordinate, regrid/remap fields. More...
 
subroutine, public register_diags_offline_transport (Time, diag, CS)
 Initialize additional diagnostics required for offline tracer transport. More...
 
subroutine, public post_offline_convergence_diags (CS, h_off, h_end, uhtr, vhtr)
 Posts diagnostics related to offline convergence diagnostics. More...
 
subroutine, public extract_offline_main (CS, uhtr, vhtr, eatr, ebtr, h_end, accumulated_time, vertical_time, dt_offline, dt_offline_vertical, skip_diffusion)
 Extracts members of the offline main control structure. All arguments are optional except the control structure itself. More...
 
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 except for the CS itself. More...
 
subroutine, public offline_transport_init (param_file, CS, diabatic_CSp, G, GV, US)
 Initializes the control structure for offline transport and reads in some of the. More...
 
subroutine read_all_input (CS)
 Coordinates the allocation and reading in all time levels of uh, vh, hend, temp, and salt from files. Used when read_all_ts_uvh. More...
 
subroutine, public offline_transport_end (CS)
 Deallocates (if necessary) arrays within the offline control structure. More...
 

Function/Subroutine Documentation

◆ extract_offline_main()

subroutine, public mom_offline_main::extract_offline_main ( type(offline_transport_cs), intent(in), target  CS,
real, dimension(:,:,:), optional, pointer  uhtr,
real, dimension(:,:,:), optional, pointer  vhtr,
real, dimension(:,:,:), optional, pointer  eatr,
real, dimension(:,:,:), optional, pointer  ebtr,
real, dimension(:,:,:), optional, pointer  h_end,
type(time_type), optional, pointer  accumulated_time,
type(time_type), optional, pointer  vertical_time,
real, intent(out), optional  dt_offline,
real, intent(out), optional  dt_offline_vertical,
logical, intent(out), optional  skip_diffusion 
)

Extracts members of the offline main control structure. All arguments are optional except the control structure itself.

Parameters
[in]csOffline control structure
uhtrRemaining zonal mass transport [H m2 ~> m3 or kg]
vhtrRemaining meridional mass transport [H m2 ~> m3 or kg]
eatrAmount of fluid entrained from the layer above within one time step [H ~> m or kg m-2]
ebtrAmount of fluid entrained from the layer below within one time step [H ~> m or kg m-2]
h_endThicknesses at the end of offline timestep [H ~> m or kg m-2]
accumulated_timeLength of time accumulated in the current offline interval
vertical_timeThe next value of accumulate_time at which to vertical processes
[out]dt_offlineTimestep used for offline tracers [T ~> s]
[out]dt_offline_verticalTimestep used for calls to tracer vertical physics [T ~> s]
[out]skip_diffusionSkips horizontal diffusion of tracers

Definition at line 1205 of file MOM_offline_main.F90.

1205  type(offline_transport_CS), target, intent(in ) :: CS !< Offline control structure
1206  ! Returned optional arguments
1207  real, dimension(:,:,:), optional, pointer :: uhtr !< Remaining zonal mass transport [H m2 ~> m3 or kg]
1208  real, dimension(:,:,:), optional, pointer :: vhtr !< Remaining meridional mass transport [H m2 ~> m3 or kg]
1209  real, dimension(:,:,:), optional, pointer :: eatr !< Amount of fluid entrained from the layer above within
1210  !! one time step [H ~> m or kg m-2]
1211  real, dimension(:,:,:), optional, pointer :: ebtr !< Amount of fluid entrained from the layer below within
1212  !! one time step [H ~> m or kg m-2]
1213  real, dimension(:,:,:), optional, pointer :: h_end !< Thicknesses at the end of offline timestep
1214  !! [H ~> m or kg m-2]
1215  type(time_type), optional, pointer :: accumulated_time !< Length of time accumulated in the
1216  !! current offline interval
1217  type(time_type), optional, pointer :: vertical_time !< The next value of accumulate_time at which to
1218  !! vertical processes
1219  real, optional, intent( out) :: dt_offline !< Timestep used for offline tracers [T ~> s]
1220  real, optional, intent( out) :: dt_offline_vertical !< Timestep used for calls to tracer
1221  !! vertical physics [T ~> s]
1222  logical, optional, intent( out) :: skip_diffusion !< Skips horizontal diffusion of tracers
1223 
1224  ! Pointers to 3d members
1225  if (present(uhtr)) uhtr => cs%uhtr
1226  if (present(vhtr)) vhtr => cs%vhtr
1227  if (present(eatr)) eatr => cs%eatr
1228  if (present(ebtr)) ebtr => cs%ebtr
1229  if (present(h_end)) h_end => cs%h_end
1230 
1231  ! Pointers to integer members which need to be modified
1232  if (present(accumulated_time)) accumulated_time => cs%accumulated_time
1233  if (present(vertical_time)) vertical_time => cs%vertical_time
1234 
1235  ! Return value of non-modified integers
1236  if (present(dt_offline)) dt_offline = cs%dt_offline
1237  if (present(dt_offline_vertical)) dt_offline_vertical = cs%dt_offline_vertical
1238  if (present(skip_diffusion)) skip_diffusion = cs%skip_diffusion
1239 

◆ insert_offline_main()

subroutine, public mom_offline_main::insert_offline_main ( type(offline_transport_cs), intent(inout)  CS,
type(ale_cs), intent(in), optional, target  ALE_CSp,
type(diabatic_cs), intent(in), optional, target  diabatic_CSp,
type(diag_ctrl), intent(in), optional, target  diag,
type(ocean_obc_type), intent(in), optional, target  OBC,
type(tracer_advect_cs), intent(in), optional, target  tracer_adv_CSp,
type(tracer_flow_control_cs), intent(in), optional, target  tracer_flow_CSp,
type(tracer_registry_type), intent(in), optional, target  tracer_Reg,
type(thermo_var_ptrs), intent(in), optional, target  tv,
type(ocean_grid_type), intent(in), optional, target  G,
type(verticalgrid_type), intent(in), optional, target  GV,
logical, intent(in), optional  x_before_y,
logical, intent(in), optional  debug 
)

Inserts (assigns values to) members of the offline main control structure. All arguments are optional except for the CS itself.

Parameters
[in,out]csOffline control structure
[in]ale_cspA pointer to the ALE control structure
[in]diabatic_cspA pointer to the diabatic control structure
[in]diagA pointer to the structure that regulates diagnostic output
[in]obcA pointer to the open boundary condition control structure
[in]tracer_adv_cspA pointer to the tracer advection control structure
[in]tracer_flow_cspA pointer to the tracer flow control control structure
[in]tracer_regA pointer to the tracer registry
[in]tvA structure pointing to various thermodynamic variables
[in]gocean grid structure
[in]gvocean vertical grid structure
[in]x_before_yIndicates which horizontal direction is advected first
[in]debugIf true, write verbose debugging messages

Definition at line 1246 of file MOM_offline_main.F90.

1246  type(offline_transport_CS), intent(inout) :: CS !< Offline control structure
1247  ! Inserted optional arguments
1248  type(ALE_CS), &
1249  target, optional, intent(in ) :: ALE_CSp !< A pointer to the ALE control structure
1250  type(diabatic_CS), &
1251  target, optional, intent(in ) :: diabatic_CSp !< A pointer to the diabatic control structure
1252  type(diag_ctrl), &
1253  target, optional, intent(in ) :: diag !< A pointer to the structure that regulates diagnostic output
1254  type(ocean_OBC_type), &
1255  target, optional, intent(in ) :: OBC !< A pointer to the open boundary condition control structure
1256  type(tracer_advect_CS), &
1257  target, optional, intent(in ) :: tracer_adv_CSp !< A pointer to the tracer advection control structure
1258  type(tracer_flow_control_CS), &
1259  target, optional, intent(in ) :: tracer_flow_CSp !< A pointer to the tracer flow control control structure
1260  type(tracer_registry_type), &
1261  target, optional, intent(in ) :: tracer_Reg !< A pointer to the tracer registry
1262  type(thermo_var_ptrs), &
1263  target, optional, intent(in ) :: tv !< A structure pointing to various thermodynamic variables
1264  type(ocean_grid_type), &
1265  target, optional, intent(in ) :: G !< ocean grid structure
1266  type(verticalGrid_type), &
1267  target, optional, intent(in ) :: GV !< ocean vertical grid structure
1268  logical, optional, intent(in ) :: x_before_y !< Indicates which horizontal direction is advected first
1269  logical, optional, intent(in ) :: debug !< If true, write verbose debugging messages
1270 
1271 
1272  if (present(ale_csp)) cs%ALE_CSp => ale_csp
1273  if (present(diabatic_csp)) cs%diabatic_CSp => diabatic_csp
1274  if (present(diag)) cs%diag => diag
1275  if (present(obc)) cs%OBC => obc
1276  if (present(tracer_adv_csp)) cs%tracer_adv_CSp => tracer_adv_csp
1277  if (present(tracer_flow_csp)) cs%tracer_flow_CSp => tracer_flow_csp
1278  if (present(tracer_reg)) cs%tracer_Reg => tracer_reg
1279  if (present(tv)) cs%tv => tv
1280  if (present(g)) cs%G => g
1281  if (present(gv)) cs%GV => gv
1282  if (present(x_before_y)) cs%x_before_y = x_before_y
1283  if (present(debug)) cs%debug = debug
1284 

◆ offline_advection_ale()

subroutine, public mom_offline_main::offline_advection_ale ( type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  Time_start,
real, intent(in)  time_interval,
type(offline_transport_cs), pointer  CS,
integer, intent(in)  id_clock_ale,
real, dimension( cs%g %isd: cs%g %ied, cs%g %jsd: cs%g %jed, cs%g %ke), intent(inout)  h_pre,
real, dimension( cs%g %isdb: cs%g %iedb, cs%g %jsd: cs%g %jed, cs%g %ke), intent(inout)  uhtr,
real, dimension( cs%g %isd: cs%g %ied, cs%g %jsdb: cs%g %jedb, cs%g %ke), intent(inout)  vhtr,
logical, intent(out)  converged 
)

3D advection is done by doing flux-limited nonlinear horizontal advection interspersed with an ALE regridding/remapping step. The loop in this routine is exited if remaining residual transports are below a runtime-specified value or a maximum number of iterations is reached.

Parameters
[in,out]fluxespointers to forcing fields
[in]time_startstarting time of a segment, as a time type
[in]time_intervaltime interval
cscontrol structure for offline module
[in]id_clock_aleClock for ALE routines
[in,out]h_prelayer thicknesses before advection
[in,out]uhtrZonal mass transport [H m2 ~> m3 or kg]
[in,out]vhtrMeridional mass transport [H m2 ~> m3 or kg]
[out]convergedTrue if the iterations have converged

Definition at line 210 of file MOM_offline_main.F90.

210  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
211  type(time_type), intent(in) :: Time_start !< starting time of a segment, as a time type
212  real, intent(in) :: time_interval !< time interval
213  type(offline_transport_CS), pointer :: CS !< control structure for offline module
214  integer, intent(in) :: id_clock_ALE !< Clock for ALE routines
215  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
216  intent(inout) :: h_pre !< layer thicknesses before advection
217  !! [H ~> m or kg m-2]
218  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
219  intent(inout) :: uhtr !< Zonal mass transport [H m2 ~> m3 or kg]
220  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), &
221  intent(inout) :: vhtr !< Meridional mass transport [H m2 ~> m3 or kg]
222  logical, intent( out) :: converged !< True if the iterations have converged
223 
224  ! Local pointers
225  type(ocean_grid_type), pointer :: G => null() ! Pointer to a structure containing
226  ! metrics and related information
227  type(verticalGrid_type), pointer :: GV => null() ! Pointer to structure containing information
228  ! about the vertical grid
229  ! Work arrays for mass transports
230  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhtr_sub
231  ! Meridional mass transports
232  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhtr_sub
233 
234  real :: prev_tot_residual, tot_residual ! Used to keep track of how close to convergence we are
235 
236  ! Variables used to keep track of layer thicknesses at various points in the code
237  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
238  h_new, &
239  h_vol
240  ! Fields for eta_diff diagnostic
241  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_pre, eta_end
242  integer :: niter, iter
243  real :: Inum_iter
244  character(len=256) :: mesg ! The text of an error message
245  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz
246  integer :: isv, iev, jsv, jev ! The valid range of the indices.
247  integer :: IsdB, IedB, JsdB, JedB
248  logical :: z_first, x_before_y
249  real :: evap_CFL_limit ! Limit on the fraction of the water that can be fluxed out of the
250  ! top layer in a timestep [nondim]
251  real :: minimum_forcing_depth ! The smallest depth over which fluxes can be applied [H ~> m or kg m-2]
252  real :: dt_iter ! The timestep to use for each iteration [T ~> s]
253 
254  integer :: nstocks
255  real :: stock_values(MAX_FIELDS_)
256  character(len=20) :: debug_msg
257  call cpu_clock_begin(cs%id_clock_offline_adv)
258 
259  ! Grid-related pointer assignments
260  g => cs%G
261  gv => cs%GV
262 
263  x_before_y = cs%x_before_y
264 
265  ! Initialize some shorthand variables from other structures
266  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
267  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
268  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
269 
270  evap_cfl_limit = cs%evap_CFL_limit
271  minimum_forcing_depth = cs%minimum_forcing_depth
272 
273  niter = cs%num_off_iter
274  inum_iter = 1./real(niter)
275  dt_iter = cs%dt_offline*inum_iter
276 
277  ! Initialize working arrays
278  h_new(:,:,:) = 0.0
279  h_vol(:,:,:) = 0.0
280  uhtr_sub(:,:,:) = 0.0
281  vhtr_sub(:,:,:) = 0.0
282 
283  ! converged should only be true if there are no remaining mass fluxes
284  converged = .false.
285 
286  ! Tracers are transported using the stored mass fluxes. Where possible, operators are Strang-split around
287  ! the call to
288  ! 1) Using the layer thicknesses and tracer concentrations from the previous timestep,
289  ! half of the accumulated vertical mixing (eatr and ebtr) is applied in the call to tracer_column_fns.
290  ! For tracers whose source/sink terms need dt, this value is set to 1/2 dt_offline
291  ! 2) Half of the accumulated surface freshwater fluxes are applied
292  !! START ITERATION
293  ! 3) Accumulated mass fluxes are used to do horizontal transport. The number of iterations used in
294  ! advect_tracer is limited to 2 (e.g x->y->x->y). The remaining mass fluxes are stored for later use
295  ! and resulting layer thicknesses fed into the next step
296  ! 4) Tracers and the h-grid are regridded and remapped in a call to ALE. This allows for layers which might
297  ! 'vanish' because of horizontal mass transport to be 'reinflated'
298  ! 5) Check that transport is done if the remaining mass fluxes equals 0 or if the max number of iterations
299  ! has been reached
300  !! END ITERATION
301  ! 6) Repeat steps 1 and 2
302  ! 7) Force a remapping to the stored layer thicknesses that correspond to the snapshot of the online model
303  ! 8) Reset T/S and h to their stored snapshotted values to prevent model drift
304 
305  ! Copy over the horizontal mass fluxes from the total mass fluxes
306  do k=1,nz ; do j=jsd,jed ; do i=isdb,iedb
307  uhtr_sub(i,j,k) = uhtr(i,j,k)
308  enddo ; enddo ; enddo
309  do k=1,nz ; do j=jsdb,jedb ; do i=isd,ied
310  vhtr_sub(i,j,k) = vhtr(i,j,k)
311  enddo ; enddo ; enddo
312  do k=1,nz ; do j=js,je ; do i=is,ie
313  h_new(i,j,k) = h_pre(i,j,k)
314  enddo ; enddo ; enddo
315 
316  if (cs%debug) then
317  call hchksum(h_pre,"h_pre before transport",g%HI)
318  call uvchksum("[uv]htr_sub before transport", uhtr_sub, vhtr_sub, g%HI)
319  endif
320  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
321  if (cs%print_adv_offline) then
322  write(mesg,'(A,ES24.16)') "Main advection starting transport: ", tot_residual
323  call mom_mesg(mesg)
324  endif
325 
326  ! This loop does essentially a flux-limited, nonlinear advection scheme until all mass fluxes
327  ! are used. ALE is done after the horizontal advection.
328  do iter=1,cs%num_off_iter
329 
330  do k=1,nz ; do j=js,je ; do i=is,ie
331  h_vol(i,j,k) = h_new(i,j,k) * g%US%L_to_m**2*g%areaT(i,j)
332  h_pre(i,j,k) = h_new(i,j,k)
333  enddo ; enddo ; enddo
334 
335  if (cs%debug) then
336  call hchksum(h_vol,"h_vol before advect",g%HI)
337  call uvchksum("[uv]htr_sub before advect", uhtr_sub, vhtr_sub, g%HI)
338  write(debug_msg, '(A,I4.4)') 'Before advect ', iter
339  call mom_tracer_chkinv(debug_msg, g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
340  endif
341 
342  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, cs%dt_offline, g, gv, cs%US, &
343  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=1, &
344  uhr_out=uhtr, vhr_out=vhtr, h_out=h_new, x_first_in=x_before_y)
345 
346  ! Switch the direction every iteration
347  x_before_y = .not. x_before_y
348 
349  ! Update the new layer thicknesses after one round of advection has happened
350  do k=1,nz ; do j=js,je ; do i=is,ie
351  h_new(i,j,k) = h_new(i,j,k) / (g%US%L_to_m**2*g%areaT(i,j))
352  enddo ; enddo ; enddo
353 
354  if (modulo(iter,cs%off_ale_mod)==0) then
355  ! Do ALE remapping/regridding to allow for more advection to occur in the next iteration
356  call pass_var(h_new,g%Domain)
357  if (cs%debug) then
358  call hchksum(h_new,"h_new before ALE",g%HI)
359  write(debug_msg, '(A,I4.4)') 'Before ALE ', iter
360  call mom_tracer_chkinv(debug_msg, g, h_new, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
361  endif
362  call cpu_clock_begin(id_clock_ale)
363  call ale_main_offline(g, gv, h_new, cs%tv, cs%tracer_Reg, cs%ALE_CSp, cs%OBC, cs%dt_offline)
364  call cpu_clock_end(id_clock_ale)
365 
366  if (cs%debug) then
367  call hchksum(h_new,"h_new after ALE",g%HI)
368  write(debug_msg, '(A,I4.4)') 'After ALE ', iter
369  call mom_tracer_chkinv(debug_msg, g, h_new, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
370  endif
371  endif
372 
373  do k=1,nz; do j=js,je ; do i=is,ie
374  uhtr_sub(i,j,k) = uhtr(i,j,k)
375  vhtr_sub(i,j,k) = vhtr(i,j,k)
376  enddo ; enddo ; enddo
377  call pass_var(h_new, g%Domain)
378  call pass_vector(uhtr_sub,vhtr_sub,g%Domain)
379 
380  ! Check for whether we've used up all the advection, or if we need to move on because
381  ! advection has stalled
382  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
383  if (cs%print_adv_offline) then
384  write(mesg,'(A,ES24.16)') "Main advection remaining transport: ", tot_residual
385  call mom_mesg(mesg)
386  endif
387  ! If all the mass transports have been used u, then quit
388  if (tot_residual == 0.0) then
389  write(mesg,*) "Converged after iteration ", iter
390  call mom_mesg(mesg)
391  converged = .true.
392  exit
393  endif
394  ! If advection has stalled or the remaining residual is less than a specified amount, quit
395  if ( (tot_residual == prev_tot_residual) .or. (tot_residual<cs%min_residual) ) then
396  converged = .false.
397  exit
398  endif
399 
400  prev_tot_residual = tot_residual
401 
402  enddo
403 
404  ! Make sure that uhtr and vhtr halos are updated
405  h_pre(:,:,:) = h_new(:,:,:)
406  call pass_vector(uhtr,vhtr,g%Domain)
407 
408  if (cs%debug) then
409  call hchksum(h_pre,"h after offline_advection_ale",g%HI)
410  call uvchksum("[uv]htr after offline_advection_ale", uhtr, vhtr, g%HI)
411  call mom_tracer_chkinv("After offline_advection_ale", g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
412  endif
413 
414  call cpu_clock_end(cs%id_clock_offline_adv)
415 

◆ offline_advection_layer()

subroutine, public mom_offline_main::offline_advection_layer ( type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  Time_start,
real, intent(in)  time_interval,
type(offline_transport_cs), pointer  CS,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  h_pre,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  eatr,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  ebtr,
real, dimension(szib_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  uhtr,
real, dimension(szi_(cs%g),szjb_(cs%g),szk_(cs%g)), intent(inout)  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.

Parameters
[in,out]fluxespointers to forcing fields
[in]time_startstarting time of a segment, as a time type
[in]time_intervalOffline transport time interval
csControl structure for offline module
[in,out]h_prelayer thicknesses before advection
[in,out]eatrEntrainment from layer above
[in,out]ebtrEntrainment from layer below
[in,out]uhtrZonal mass transport
[in,out]vhtrMeridional mass transport

Definition at line 846 of file MOM_offline_main.F90.

846  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
847  type(time_type), intent(in) :: Time_start !< starting time of a segment, as a time type
848  real, intent(in) :: time_interval !< Offline transport time interval
849  type(offline_transport_CS), pointer :: CS !< Control structure for offline module
850  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_pre !< layer thicknesses before advection
851  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: eatr !< Entrainment from layer above
852  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: ebtr !< Entrainment from layer below
853  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: uhtr !< Zonal mass transport
854  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(inout) :: vhtr !< Meridional mass transport
855  ! Local pointers
856  type(ocean_grid_type), pointer :: G => null() ! Pointer to a structure containing
857  ! metrics and related information
858  type(verticalGrid_type), pointer :: GV => null() ! Pointer to structure containing information
859  ! about the vertical grid
860  ! Remaining zonal mass transports
861  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhtr_sub
862  ! Remaining meridional mass transports
863  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhtr_sub
864 
865  real :: sum_abs_fluxes, sum_u, sum_v ! Used to keep track of how close to convergence we are
866  real :: dt_offline
867 
868  ! Local variables
869  ! Vertical diffusion related variables
870  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
871  eatr_sub, &
872  ebtr_sub
873  ! Variables used to keep track of layer thicknesses at various points in the code
874  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
875  h_new, &
876  h_vol
877  ! Work arrays for temperature and salinity
878  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
879  temp_old, salt_old, &
880  temp_mean, salt_mean, &
881  zero_3dh !
882  integer :: niter, iter
883  real :: Inum_iter
884  real :: dt_iter ! The timestep of each iteration [T ~> s]
885  logical :: converged
886  character(len=160) :: mesg ! The text of an error message
887  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz
888  integer :: isv, iev, jsv, jev ! The valid range of the indices.
889  integer :: IsdB, IedB, JsdB, JedB
890  logical :: z_first, x_before_y
891 
892  g => cs%G ; gv => cs%GV
893  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
894  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
895  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
896 
897  dt_iter = cs%US%s_to_T * time_interval / real(max(1, cs%num_off_iter))
898 
899  do iter=1,cs%num_off_iter
900 
901  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
902  eatr_sub(i,j,k) = eatr(i,j,k)
903  ebtr_sub(i,j,k) = ebtr(i,j,k)
904  enddo ; enddo ; enddo
905 
906  do k = 1, nz ; do j=js-1,je+1 ; do i=is-2,ie+1
907  uhtr_sub(i,j,k) = uhtr(i,j,k)
908  enddo ; enddo ; enddo
909 
910  do k = 1, nz ; do j=js-2,je+1 ; do i=is-1,ie+1
911  vhtr_sub(i,j,k) = vhtr(i,j,k)
912  enddo ; enddo ; enddo
913 
914 
915  ! Calculate 3d mass transports to be used in this iteration
916  call limit_mass_flux_3d(g, gv, uhtr_sub, vhtr_sub, eatr_sub, ebtr_sub, h_pre)
917 
918  if (z_first) then
919  ! First do vertical advection
920  call update_h_vertical_flux(g, gv, eatr_sub, ebtr_sub, h_pre, h_new)
921  call call_tracer_column_fns(h_pre, h_new, eatr_sub, ebtr_sub, &
922  fluxes, cs%mld, dt_iter, g, gv, cs%US, cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
923  ! We are now done with the vertical mass transports, so now h_new is h_sub
924  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
925  h_pre(i,j,k) = h_new(i,j,k)
926  enddo ; enddo ; enddo
927  call pass_var(h_pre,g%Domain)
928 
929  ! Second zonal and meridional advection
930  call update_h_horizontal_flux(g, gv, uhtr_sub, vhtr_sub, h_pre, h_new)
931  do k = 1, nz ; do i = is-1, ie+1 ; do j=js-1, je+1
932  h_vol(i,j,k) = h_pre(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)
933  enddo ; enddo ; enddo
934  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, dt_iter, g, gv, cs%US, &
935  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=30, x_first_in=x_before_y)
936 
937  ! Done with horizontal so now h_pre should be h_new
938  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
939  h_pre(i,j,k) = h_new(i,j,k)
940  enddo ; enddo ; enddo
941 
942  endif
943 
944  if (.not. z_first) then
945 
946  ! First zonal and meridional advection
947  call update_h_horizontal_flux(g, gv, uhtr_sub, vhtr_sub, h_pre, h_new)
948  do k = 1, nz ; do i = is-1, ie+1 ; do j=js-1, je+1
949  h_vol(i,j,k) = h_pre(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)
950  enddo ; enddo ; enddo
951  call advect_tracer(h_pre, uhtr_sub, vhtr_sub, cs%OBC, dt_iter, g, gv, cs%US, &
952  cs%tracer_adv_CSp, cs%tracer_Reg, h_vol, max_iter_in=30, x_first_in=x_before_y)
953 
954  ! Done with horizontal so now h_pre should be h_new
955  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
956  h_pre(i,j,k) = h_new(i,j,k)
957  enddo ; enddo ; enddo
958 
959  ! Second vertical advection
960  call update_h_vertical_flux(g, gv, eatr_sub, ebtr_sub, h_pre, h_new)
961  call call_tracer_column_fns(h_pre, h_new, eatr_sub, ebtr_sub, &
962  fluxes, cs%mld, dt_iter, g, gv, cs%US, cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
963  ! We are now done with the vertical mass transports, so now h_new is h_sub
964  do k = 1, nz ; do i=is-1,ie+1 ; do j=js-1,je+1
965  h_pre(i,j,k) = h_new(i,j,k)
966  enddo ; enddo ; enddo
967 
968  endif
969 
970  ! Update remaining transports
971  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
972  eatr(i,j,k) = eatr(i,j,k) - eatr_sub(i,j,k)
973  ebtr(i,j,k) = ebtr(i,j,k) - ebtr_sub(i,j,k)
974  enddo ; enddo ; enddo
975 
976  do k = 1, nz ; do j=js-1,je+1 ; do i=is-2,ie+1
977  uhtr(i,j,k) = uhtr(i,j,k) - uhtr_sub(i,j,k)
978  enddo ; enddo ; enddo
979 
980  do k = 1, nz ; do j=js-2,je+1 ; do i=is-1,ie+1
981  vhtr(i,j,k) = vhtr(i,j,k) - vhtr_sub(i,j,k)
982  enddo ; enddo ; enddo
983 
984  call pass_var(eatr,g%Domain)
985  call pass_var(ebtr,g%Domain)
986  call pass_var(h_pre,g%Domain)
987  call pass_vector(uhtr,vhtr,g%Domain)
988  !
989  ! Calculate how close we are to converging by summing the remaining fluxes at each point
990  sum_abs_fluxes = 0.0
991  sum_u = 0.0
992  sum_v = 0.0
993  do k=1,nz; do j=js,je; do i=is,ie
994  sum_u = sum_u + abs(uhtr(i-1,j,k))+abs(uhtr(i,j,k))
995  sum_v = sum_v + abs(vhtr(i,j-1,k))+abs(vhtr(i,j,k))
996  sum_abs_fluxes = sum_abs_fluxes + abs(eatr(i,j,k)) + abs(ebtr(i,j,k)) + abs(uhtr(i-1,j,k)) + &
997  abs(uhtr(i,j,k)) + abs(vhtr(i,j-1,k)) + abs(vhtr(i,j,k))
998  enddo ; enddo ; enddo
999  call sum_across_pes(sum_abs_fluxes)
1000 
1001  write(mesg,*) "offline_advection_layer: Remaining u-flux, v-flux:", sum_u, sum_v
1002  call mom_mesg(mesg)
1003  if (sum_abs_fluxes==0) then
1004  write(mesg,*) 'offline_advection_layer: Converged after iteration', iter
1005  call mom_mesg(mesg)
1006  exit
1007  endif
1008 
1009  ! Switch order of Strang split every iteration
1010  z_first = .not. z_first
1011  x_before_y = .not. x_before_y
1012  enddo
1013 

◆ offline_diabatic_ale()

subroutine, public mom_offline_main::offline_diabatic_ale ( type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  Time_start,
type(time_type), intent(in)  Time_end,
type(offline_transport_cs), pointer  CS,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  h_pre,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  eatr,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  ebtr 
)

The vertical/diabatic driver for offline tracers. First the eatr/ebtr associated with the interpolated vertical diffusivities are calculated and then any tracer column functions are done which can include vertical diffuvities and source/sink terms.

Parameters
[in,out]fluxespointers to forcing fields
[in]time_startstarting time of a segment, as a time type
[in]time_endtime interval
cscontrol structure from initialize_MOM
[in,out]h_prelayer thicknesses before advection [H ~> m or kg m-2]
[in,out]eatrEntrainment from layer above [H ~> m or kg m-2]
[in,out]ebtrEntrainment from layer below [H ~> m or kg m-2]

Definition at line 654 of file MOM_offline_main.F90.

654 
655  type(forcing), intent(inout) :: fluxes !< pointers to forcing fields
656  type(time_type), intent(in) :: Time_start !< starting time of a segment, as a time type
657  type(time_type), intent(in) :: Time_end !< time interval
658  type(offline_transport_CS), pointer :: CS !< control structure from initialize_MOM
659  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
660  intent(inout) :: h_pre !< layer thicknesses before advection [H ~> m or kg m-2]
661  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
662  intent(inout) :: eatr !< Entrainment from layer above [H ~> m or kg m-2]
663  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
664  intent(inout) :: ebtr !< Entrainment from layer below [H ~> m or kg m-2]
665 
666  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: &
667  sw, sw_vis, sw_nir !< Save old values of shortwave radiation [Q R Z T-1 ~> W m-2]
668  real :: hval
669  integer :: i,j,k
670  integer :: is, ie, js, je, nz
671  integer :: k_nonzero
672  real :: stock_values(MAX_FIELDS_)
673  real :: Kd_bot
674  integer :: nstocks
675  nz = cs%GV%ke
676  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec
677 
678  call cpu_clock_begin(cs%id_clock_offline_diabatic)
679 
680  call mom_mesg("Applying tracer source, sinks, and vertical mixing")
681 
682  if (cs%debug) then
683  call hchksum(h_pre,"h_pre before offline_diabatic_ale",cs%G%HI)
684  call hchksum(eatr,"eatr before offline_diabatic_ale",cs%G%HI)
685  call hchksum(ebtr,"ebtr before offline_diabatic_ale",cs%G%HI)
686  call mom_tracer_chkinv("Before offline_diabatic_ale", cs%G, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
687  endif
688 
689  eatr(:,:,:) = 0.
690  ebtr(:,:,:) = 0.
691  ! Calculate eatr and ebtr if vertical diffusivity is read
692  ! Because the saved remapped diagnostics from the online run assume a zero minimum thickness
693  ! but ALE may have a minimum thickness. Flood the diffusivities for all layers with the value
694  ! of Kd closest to the bottom which is non-zero
695  do j=js,je ; do i=is,ie
696  k_nonzero = nz+1
697  ! Find the nonzero bottom Kd
698  do k=nz+1,1,-1
699  if (cs%Kd(i,j,k)>0.) then
700  kd_bot = cs%Kd(i,j,k)
701  k_nonzero = k
702  exit
703  endif
704  enddo
705  ! Flood the bottom interfaces
706  do k=k_nonzero,nz+1
707  cs%Kd(i,j,k) = kd_bot
708  enddo
709  enddo ; enddo
710 
711  do j=js,je ; do i=is,ie
712  eatr(i,j,1) = 0.
713  enddo ; enddo
714  do k=2,nz ; do j=js,je ; do i=is,ie
715  hval=1.0/(cs%GV%H_subroundoff + 0.5*(h_pre(i,j,k-1) + h_pre(i,j,k)))
716  eatr(i,j,k) = (cs%GV%m_to_H**2*cs%US%T_to_s) * cs%dt_offline_vertical * hval * cs%Kd(i,j,k)
717  ebtr(i,j,k-1) = eatr(i,j,k)
718  enddo ; enddo ; enddo
719  do j=js,je ; do i=is,ie
720  ebtr(i,j,nz) = 0.
721  enddo ; enddo
722 
723  ! Add diurnal cycle for shortwave radiation (only used if run in ocean-only mode)
724  if (cs%diurnal_SW .and. cs%read_sw) then
725  sw(:,:) = fluxes%sw(:,:)
726  sw_vis(:,:) = fluxes%sw_vis_dir(:,:)
727  sw_nir(:,:) = fluxes%sw_nir_dir(:,:)
728  call offline_add_diurnal_sw(fluxes, cs%G, time_start, time_end)
729  endif
730 
731  if (associated(cs%optics)) &
732  call set_pen_shortwave(cs%optics, fluxes, cs%G, cs%GV, cs%US, cs%diabatic_aux_CSp, &
733  cs%opacity_CSp, cs%tracer_flow_CSp)
734 
735  ! Note that tracerBoundaryFluxesInOut within this subroutine should NOT be called
736  ! as the freshwater fluxes have already been accounted for
737  call call_tracer_column_fns(h_pre, h_pre, eatr, ebtr, fluxes, cs%MLD, cs%dt_offline_vertical, &
738  cs%G, cs%GV, cs%US, cs%tv, cs%optics, cs%tracer_flow_CSp, cs%debug)
739 
740  if (cs%diurnal_SW .and. cs%read_sw) then
741  fluxes%sw(:,:) = sw(:,:)
742  fluxes%sw_vis_dir(:,:) = sw_vis(:,:)
743  fluxes%sw_nir_dir(:,:) = sw_nir(:,:)
744  endif
745 
746  if (cs%debug) then
747  call hchksum(h_pre,"h_pre after offline_diabatic_ale",cs%G%HI)
748  call hchksum(eatr,"eatr after offline_diabatic_ale",cs%G%HI)
749  call hchksum(ebtr,"ebtr after offline_diabatic_ale",cs%G%HI)
750  call mom_tracer_chkinv("After offline_diabatic_ale", cs%G, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
751  endif
752 
753  call cpu_clock_end(cs%id_clock_offline_diabatic)
754 

◆ offline_fw_fluxes_into_ocean()

subroutine, public mom_offline_main::offline_fw_fluxes_into_ocean ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(offline_transport_cs), intent(inout)  CS,
type(forcing), intent(inout)  fluxes,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(inout)  h,
real, dimension(szi_(g),szj_(g)), intent(in), optional  in_flux_optional 
)

Apply positive freshwater fluxes (into the ocean) and update netMassOut with only the negative (out of the ocean) fluxes.

Parameters
[in,out]csOffline control structure
[in]gGrid structure
[in]gvocean vertical grid structure
[in,out]fluxesSurface fluxes container
[in,out]hLayer thickness [H ~> m or kg m-2]
[in]in_flux_optionalThe total time-integrated amount

Definition at line 760 of file MOM_offline_main.F90.

760  type(offline_transport_CS), intent(inout) :: CS !< Offline control structure
761  type(ocean_grid_type), intent(in) :: G !< Grid structure
762  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
763  type(forcing), intent(inout) :: fluxes !< Surface fluxes container
764  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
765  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
766  real, dimension(SZI_(G),SZJ_(G)), &
767  optional, intent(in) :: in_flux_optional !< The total time-integrated amount
768  !! of tracer that leaves with freshwater
769 
770  integer :: i, j, m
771  real, dimension(SZI_(G),SZJ_(G)) :: negative_fw !< store all negative fluxes
772  logical :: update_h !< Flag for whether h should be updated
773 
774  if ( present(in_flux_optional) ) &
775  call mom_error(warning, "Positive freshwater fluxes with non-zero tracer concentration not supported yet")
776 
777  ! Set all fluxes to 0
778  negative_fw(:,:) = 0.
779 
780  ! Sort fluxes into positive and negative
781  do j=g%jsc,g%jec ; do i=g%isc,g%iec
782  if (fluxes%netMassOut(i,j)<0.0) then
783  negative_fw(i,j) = fluxes%netMassOut(i,j)
784  fluxes%netMassOut(i,j) = 0.
785  endif
786  enddo ; enddo
787 
788  if (cs%debug) then
789  call hchksum(h,"h before fluxes into ocean",g%HI)
790  call mom_tracer_chkinv("Before fluxes into ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
791  endif
792  do m = 1,cs%tracer_reg%ntr
793  ! Layer thicknesses should only be updated after the last tracer is finished
794  update_h = ( m == cs%tracer_reg%ntr )
795  call applytracerboundaryfluxesinout(g, gv, cs%tracer_reg%tr(m)%t, cs%dt_offline, fluxes, h, &
796  cs%evap_CFL_limit, cs%minimum_forcing_depth, update_h_opt = update_h)
797  enddo
798  if (cs%debug) then
799  call hchksum(h,"h after fluxes into ocean",g%HI)
800  call mom_tracer_chkinv("After fluxes into ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
801  endif
802 
803  ! Now that fluxes into the ocean are done, save the negative fluxes for later
804  fluxes%netMassOut(:,:) = negative_fw(:,:)
805 

◆ offline_fw_fluxes_out_ocean()

subroutine, public mom_offline_main::offline_fw_fluxes_out_ocean ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(offline_transport_cs), intent(inout)  CS,
type(forcing), intent(inout)  fluxes,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(inout)  h,
real, dimension(szi_(g),szj_(g)), intent(in), optional  out_flux_optional 
)

Apply negative freshwater fluxes (out of the ocean)

Parameters
[in,out]csOffline control structure
[in]gGrid structure
[in]gvocean vertical grid structure
[in,out]fluxesSurface fluxes container
[in,out]hLayer thickness [H ~> m or kg m-2]
[in]out_flux_optionalThe total time-integrated amount

Definition at line 810 of file MOM_offline_main.F90.

810  type(offline_transport_CS), intent(inout) :: CS !< Offline control structure
811  type(ocean_grid_type), intent(in) :: G !< Grid structure
812  type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
813  type(forcing), intent(inout) :: fluxes !< Surface fluxes container
814  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
815  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
816  real, dimension(SZI_(G),SZJ_(G)), &
817  optional, intent(in) :: out_flux_optional !< The total time-integrated amount
818  !! of tracer that leaves with freshwater
819 
820  integer :: m
821  logical :: update_h !< Flag for whether h should be updated
822 
823  if ( present(out_flux_optional) ) &
824  call mom_error(warning, "Negative freshwater fluxes with non-zero tracer concentration not supported yet")
825 
826  if (cs%debug) then
827  call hchksum(h,"h before fluxes out of ocean",g%HI)
828  call mom_tracer_chkinv("Before fluxes out of ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
829  endif
830  do m = 1, cs%tracer_reg%ntr
831  ! Layer thicknesses should only be updated after the last tracer is finished
832  update_h = ( m == cs%tracer_reg%ntr )
833  call applytracerboundaryfluxesinout(g, gv, cs%tracer_reg%tr(m)%t, cs%dt_offline, fluxes, h, &
834  cs%evap_CFL_limit, cs%minimum_forcing_depth, update_h_opt = update_h)
835  enddo
836  if (cs%debug) then
837  call hchksum(h,"h after fluxes out of ocean",g%HI)
838  call mom_tracer_chkinv("Before fluxes out of ocean", g, h, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
839  endif
840 

◆ offline_redistribute_residual()

subroutine, public mom_offline_main::offline_redistribute_residual ( type(offline_transport_cs), pointer  CS,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  h_pre,
real, dimension(szib_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  uhtr,
real, dimension(szi_(cs%g),szjb_(cs%g),szk_(cs%g)), intent(inout)  vhtr,
logical, intent(in)  converged 
)

In the case where the main advection routine did not converge, something needs to be done with the remaining transport. Two different ways are offered, 'barotropic' means that the residual is distributed equally throughout the water column. 'upwards' attempts to redistribute the transport in the layers above and will eventually work down the entire water column.

Parameters
cscontrol structure from initialize_MOM
[in,out]h_prelayer thicknesses before advection
[in,out]uhtrZonal mass transport
[in,out]vhtrMeridional mass transport
[in]convergedTrue if the iterations have converged

Definition at line 423 of file MOM_offline_main.F90.

423  type(offline_transport_CS), pointer :: CS !< control structure from initialize_MOM
424  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
425  intent(inout) :: h_pre !< layer thicknesses before advection
426  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), &
427  intent(inout) :: uhtr !< Zonal mass transport
428  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), &
429  intent(inout) :: vhtr !< Meridional mass transport
430  logical, intent(in ) :: converged !< True if the iterations have converged
431 
432  type(ocean_grid_type), pointer :: G => null() ! Pointer to a structure containing
433  ! metrics and related information
434  type(verticalGrid_type), pointer :: GV => null() ! Pointer to structure containing information
435  ! about the vertical grid
436  logical :: x_before_y
437  ! Variables used to keep track of layer thicknesses at various points in the code
438  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: &
439  h_new, &
440  h_vol
441 
442  ! Used to calculate the eta diagnostics
443  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_work
444  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: uhr !< Zonal mass transport
445  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)) :: vhr !< Meridional mass transport
446 
447  character(len=256) :: mesg ! The text of an error message
448  integer :: i, j, k, m, is, ie, js, je, isd, ied, jsd, jed, nz, iter
449  real :: prev_tot_residual, tot_residual, stock_values(MAX_FIELDS_)
450  integer :: nstocks
451 
452  ! Assign grid pointers
453  g => cs%G
454  gv => cs%GV
455 
456  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
457  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
458 
459  x_before_y = cs%x_before_y
460 
461  if (cs%id_eta_pre_distribute>0) then
462  eta_work(:,:) = 0.0
463  do k=1,nz ; do j=js,je ; do i=is,ie
464  if (h_pre(i,j,k)>gv%Angstrom_H) then
465  eta_work(i,j) = eta_work(i,j) + h_pre(i,j,k)
466  endif
467  enddo ; enddo ; enddo
468  call post_data(cs%id_eta_pre_distribute,eta_work,cs%diag)
469  endif
470 
471  ! These are used to find out how much will be redistributed in this routine
472  if (cs%id_h_redist>0) call post_data(cs%id_h_redist, h_pre, cs%diag)
473  if (cs%id_uhr_redist>0) call post_data(cs%id_uhr_redist, uhtr, cs%diag)
474  if (cs%id_vhr_redist>0) call post_data(cs%id_vhr_redist, vhtr, cs%diag)
475 
476  if (converged) return
477 
478  if (cs%debug) then
479  call mom_tracer_chkinv("Before redistribute ", g, h_pre, cs%tracer_reg%Tr, cs%tracer_reg%ntr)
480  endif
481 
482  call cpu_clock_begin(cs%id_clock_redistribute)
483 
484  if (cs%redistribute_upwards .or. cs%redistribute_barotropic) then
485  do iter = 1, cs%num_off_iter
486 
487  ! Perform upwards redistribution
488  if (cs%redistribute_upwards) then
489 
490  ! Calculate the layer volumes at beginning of redistribute
491  do k=1,nz ; do j=js,je ; do i=is,ie
492  h_vol(i,j,k) = h_pre(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)
493  enddo ; enddo ; enddo
494  call pass_var(h_vol,g%Domain)
495  call pass_vector(uhtr,vhtr,g%Domain)
496 
497  ! Store volumes for advect_tracer
498  h_pre(:,:,:) = h_vol(:,:,:)
499 
500  if (cs%debug) then
501  call mom_tracer_chksum("Before upwards redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
502  call uvchksum("[uv]tr before upwards redistribute", uhtr, vhtr, g%HI)
503  endif
504 
505  if (x_before_y) then
506  call distribute_residual_uh_upwards(g, gv, h_vol, uhtr)
507  call distribute_residual_vh_upwards(g, gv, h_vol, vhtr)
508  else
509  call distribute_residual_vh_upwards(g, gv, h_vol, vhtr)
510  call distribute_residual_uh_upwards(g, gv, h_vol, uhtr)
511  endif
512 
513  call advect_tracer(h_pre, uhtr, vhtr, cs%OBC, cs%dt_offline, g, gv, cs%US, &
514  cs%tracer_adv_CSp, cs%tracer_Reg, h_prev_opt = h_pre, max_iter_in=1, &
515  h_out=h_new, uhr_out=uhr, vhr_out=vhr, x_first_in=x_before_y)
516 
517  if (cs%debug) then
518  call mom_tracer_chksum("After upwards redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
519  endif
520 
521  ! Convert h_new back to layer thickness for ALE remapping
522  do k=1,nz ; do j=js,je ; do i=is,ie
523  uhtr(i,j,k) = uhr(i,j,k)
524  vhtr(i,j,k) = vhr(i,j,k)
525  h_vol(i,j,k) = h_new(i,j,k)
526  h_new(i,j,k) = h_new(i,j,k) / (g%US%L_to_m**2*g%areaT(i,j))
527  h_pre(i,j,k) = h_new(i,j,k)
528  enddo ; enddo ; enddo
529 
530  endif ! redistribute upwards
531 
532  ! Perform barotropic redistribution
533  if (cs%redistribute_barotropic) then
534 
535  ! Calculate the layer volumes at beginning of redistribute
536  do k=1,nz ; do j=js,je ; do i=is,ie
537  h_vol(i,j,k) = h_pre(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)
538  enddo ; enddo ; enddo
539  call pass_var(h_vol,g%Domain)
540  call pass_vector(uhtr,vhtr,g%Domain)
541 
542  ! Copy h_vol to h_pre for advect_tracer routine
543  h_pre(:,:,:) = h_vol(:,:,:)
544 
545  if (cs%debug) then
546  call mom_tracer_chksum("Before barotropic redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
547  call uvchksum("[uv]tr before upwards redistribute", uhtr, vhtr, g%HI)
548  endif
549 
550  if (x_before_y) then
551  call distribute_residual_uh_barotropic(g, gv, h_vol, uhtr)
552  call distribute_residual_vh_barotropic(g, gv, h_vol, vhtr)
553  else
554  call distribute_residual_vh_barotropic(g, gv, h_vol, vhtr)
555  call distribute_residual_uh_barotropic(g, gv, h_vol, uhtr)
556  endif
557 
558  call advect_tracer(h_pre, uhtr, vhtr, cs%OBC, cs%dt_offline, g, gv, cs%US, &
559  cs%tracer_adv_CSp, cs%tracer_Reg, h_prev_opt = h_pre, max_iter_in=1, &
560  h_out=h_new, uhr_out=uhr, vhr_out=vhr, x_first_in=x_before_y)
561 
562  if (cs%debug) then
563  call mom_tracer_chksum("After barotropic redistribute ", cs%tracer_Reg%Tr, cs%tracer_Reg%ntr, g)
564  endif
565 
566  ! Convert h_new back to layer thickness for ALE remapping
567  do k=1,nz ; do j=js,je ; do i=is,ie
568  uhtr(i,j,k) = uhr(i,j,k)
569  vhtr(i,j,k) = vhr(i,j,k)
570  h_vol(i,j,k) = h_new(i,j,k)
571  h_new(i,j,k) = h_new(i,j,k) / (g%US%L_to_m**2*g%areaT(i,j))
572  h_pre(i,j,k) = h_new(i,j,k)
573  enddo ; enddo ; enddo
574 
575  endif ! redistribute barotropic
576 
577  ! Check to see if all transport has been exhausted
578  tot_residual = remaining_transport_sum(cs, uhtr, vhtr)
579  if (cs%print_adv_offline) then
580  write(mesg,'(A,ES24.16)') "Residual advection remaining transport: ", tot_residual
581  call mom_mesg(mesg)
582  endif
583  ! If the remaining residual is 0, then this return is done
584  if (tot_residual==0.0 ) then
585  exit
586  endif
587 
588  if ( (tot_residual == prev_tot_residual) .or. (tot_residual<cs%min_residual) ) exit
589  prev_tot_residual = tot_residual
590 
591  enddo ! Redistribution iteration
592  endif ! If one of the redistribution routines is requested
593 
594  if (cs%id_eta_post_distribute>0) then
595  eta_work(:,:) = 0.0
596  do k=1,nz ; do j=js,je ; do i=is,ie
597  if (h_pre(i,j,k)>gv%Angstrom_H) then
598  eta_work(i,j) = eta_work(i,j) + h_pre(i,j,k)
599  endif
600  enddo ; enddo ; enddo
601  call post_data(cs%id_eta_post_distribute,eta_work,cs%diag)
602  endif
603 
604  if (cs%id_uhr>0) call post_data(cs%id_uhr,uhtr,cs%diag)
605  if (cs%id_vhr>0) call post_data(cs%id_vhr,vhtr,cs%diag)
606 
607  if (cs%debug) then
608  call hchksum(h_pre,"h_pre after redistribute",g%HI)
609  call uvchksum("uhtr after redistribute", uhtr, vhtr, g%HI)
610  call mom_tracer_chkinv("after redistribute ", g, h_new, cs%tracer_Reg%Tr, cs%tracer_Reg%ntr)
611  endif
612 
613  call cpu_clock_end(cs%id_clock_redistribute)
614 

◆ offline_transport_end()

subroutine, public mom_offline_main::offline_transport_end ( type(offline_transport_cs), pointer  CS)

Deallocates (if necessary) arrays within the offline control structure.

Parameters
csControl structure for offline module

Definition at line 1506 of file MOM_offline_main.F90.

1506  type(offline_transport_CS), pointer :: CS !< Control structure for offline module
1507 
1508  ! Explicitly allocate all allocatable arrays
1509  deallocate(cs%uhtr)
1510  deallocate(cs%vhtr)
1511  deallocate(cs%eatr)
1512  deallocate(cs%ebtr)
1513  deallocate(cs%h_end)
1514  deallocate(cs%netMassOut)
1515  deallocate(cs%netMassIn)
1516  deallocate(cs%Kd)
1517  if (cs%read_mld) deallocate(cs%mld)
1518  if (cs%read_all_ts_uvh) then
1519  deallocate(cs%uhtr_all)
1520  deallocate(cs%vhtr_all)
1521  deallocate(cs%hend_all)
1522  deallocate(cs%temp_all)
1523  deallocate(cs%salt_all)
1524  endif
1525 
1526  deallocate(cs)
1527 

◆ offline_transport_init()

subroutine, public mom_offline_main::offline_transport_init ( type(param_file_type), intent(in)  param_file,
type(offline_transport_cs), pointer  CS,
type(diabatic_cs), intent(in)  diabatic_CSp,
type(ocean_grid_type), intent(in), target  G,
type(verticalgrid_type), intent(in), target  GV,
type(unit_scale_type), intent(in), target  US 
)

Initializes the control structure for offline transport and reads in some of the.

Parameters
[in]param_fileA structure to parse for run-time parameters
csOffline control structure
[in]diabatic_cspThe diabatic control structure
[in]gocean grid structure
[in]gvocean vertical grid structure
[in]usA dimensional unit scaling type

Definition at line 1290 of file MOM_offline_main.F90.

1290 
1291  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1292  type(offline_transport_CS), pointer :: CS !< Offline control structure
1293  type(diabatic_CS), intent(in) :: diabatic_CSp !< The diabatic control structure
1294  type(ocean_grid_type), target, intent(in) :: G !< ocean grid structure
1295  type(verticalGrid_type), target, intent(in) :: GV !< ocean vertical grid structure
1296  type(unit_scale_type), target, intent(in) :: US !< A dimensional unit scaling type
1297 
1298  character(len=40) :: mdl = "offline_transport"
1299  character(len=20) :: redistribute_method
1300 
1301  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
1302  integer :: IsdB, IedB, JsdB, JedB
1303 
1304  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
1305  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1306  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1307 
1308  call calltree_enter("offline_transport_init, MOM_offline_control.F90")
1309 
1310  if (associated(cs)) then
1311  call mom_error(warning, "offline_transport_init called with an associated "// &
1312  "control structure.")
1313  return
1314  endif
1315  allocate(cs)
1316  call log_version(param_file, mdl,version, "This module allows for tracers to be run offline")
1317 
1318  ! Determining the internal unit scaling factors for this run.
1319  cs%US => us
1320 
1321  ! Parse MOM_input for offline control
1322  call get_param(param_file, mdl, "OFFLINEDIR", cs%offlinedir, &
1323  "Input directory where the offline fields can be found", fail_if_missing = .true.)
1324  call get_param(param_file, mdl, "OFF_SUM_FILE", cs%sum_file, &
1325  "Filename where the accumulated fields can be found", fail_if_missing = .true.)
1326  call get_param(param_file, mdl, "OFF_SNAP_FILE", cs%snap_file, &
1327  "Filename where snapshot fields can be found", fail_if_missing = .true.)
1328  call get_param(param_file, mdl, "OFF_MEAN_FILE", cs%mean_file, &
1329  "Filename where averaged fields can be found", fail_if_missing = .true.)
1330  call get_param(param_file, mdl, "OFF_SURF_FILE", cs%surf_file, &
1331  "Filename where averaged fields can be found", fail_if_missing = .true.)
1332  call get_param(param_file, mdl, "NUMTIME", cs%numtime, &
1333  "Number of timelevels in offline input files", fail_if_missing = .true.)
1334  call get_param(param_file, mdl, "NK_INPUT", cs%nk_input, &
1335  "Number of vertical levels in offline input files", default = nz)
1336  call get_param(param_file, mdl, "DT_OFFLINE", cs%dt_offline, &
1337  "Length of time between reading in of input fields", units='s', scale=us%s_to_T, fail_if_missing = .true.)
1338  call get_param(param_file, mdl, "DT_OFFLINE_VERTICAL", cs%dt_offline_vertical, &
1339  "Length of the offline timestep for tracer column sources/sinks " //&
1340  "This should be set to the length of the coupling timestep for " //&
1341  "tracers which need shortwave fluxes", units="s", scale=us%s_to_T, fail_if_missing = .true.)
1342  call get_param(param_file, mdl, "START_INDEX", cs%start_index, &
1343  "Which time index to start from", default=1)
1344  call get_param(param_file, mdl, "FIELDS_ARE_OFFSET", cs%fields_are_offset, &
1345  "True if the time-averaged fields and snapshot fields "//&
1346  "are offset by one time level", default=.false.)
1347  call get_param(param_file, mdl, "REDISTRIBUTE_METHOD", redistribute_method, &
1348  "Redistributes any remaining horizontal fluxes throughout " //&
1349  "the rest of water column. Options are 'barotropic' which " //&
1350  "evenly distributes flux throughout the entire water column, " //&
1351  "'upwards' which adds the maximum of the remaining flux in " //&
1352  "each layer above, both which first applies upwards and then " //&
1353  "barotropic, and 'none' which does no redistribution", &
1354  default='barotropic')
1355  call get_param(param_file, mdl, "NUM_OFF_ITER", cs%num_off_iter, &
1356  "Number of iterations to subdivide the offline tracer advection and diffusion", &
1357  default = 60)
1358  call get_param(param_file, mdl, "OFF_ALE_MOD", cs%off_ale_mod, &
1359  "Sets how many horizontal advection steps are taken before an ALE " //&
1360  "remapping step is done. 1 would be x->y->ALE, 2 would be" //&
1361  "x->y->x->y->ALE", default = 1)
1362  call get_param(param_file, mdl, "PRINT_ADV_OFFLINE", cs%print_adv_offline, &
1363  "Print diagnostic output every advection subiteration",default=.false.)
1364  call get_param(param_file, mdl, "SKIP_DIFFUSION_OFFLINE", cs%skip_diffusion, &
1365  "Do not do horizontal diffusion",default=.false.)
1366  call get_param(param_file, mdl, "READ_SW", cs%read_sw, &
1367  "Read in shortwave radiation field instead of using values from the coupler"//&
1368  "when in offline tracer mode",default=.false.)
1369  call get_param(param_file, mdl, "READ_MLD", cs%read_mld, &
1370  "Read in mixed layer depths for tracers which exchange with the atmosphere"//&
1371  "when in offline tracer mode",default=.false.)
1372  call get_param(param_file, mdl, "MLD_VAR_NAME", cs%mld_var_name, &
1373  "Name of the variable containing the depth of active mixing",&
1374  default='ePBL_h_ML')
1375  call get_param(param_file, mdl, "OFFLINE_ADD_DIURNAL_SW", cs%diurnal_sw, &
1376  "Adds a synthetic diurnal cycle in the same way that the ice " // &
1377  "model would have when time-averaged fields of shortwave " // &
1378  "radiation are read in", default=.false.)
1379  call get_param(param_file, mdl, "KD_MAX", cs%Kd_max, &
1380  "The maximum permitted increment for the diapycnal "//&
1381  "diffusivity from TKE-based parameterizations, or a "//&
1382  "negative value for no limit.", units="m2 s-1", default=-1.0)
1383  call get_param(param_file, mdl, "MIN_RESIDUAL_TRANSPORT", cs%min_residual, &
1384  "How much remaining transport before the main offline advection "// &
1385  "is exited. The default value corresponds to about 1 meter of " // &
1386  "difference in a grid cell", default = 1.e9)
1387  call get_param(param_file, mdl, "READ_ALL_TS_UVH", cs%read_all_ts_uvh, &
1388  "Reads all time levels of a subset of the fields necessary to run " // &
1389  "the model offline. This can require a large amount of memory "// &
1390  "and will make initialization very slow. However, for offline "// &
1391  "runs spanning more than a year this can reduce total I/O overhead", &
1392  default = .false.)
1393 
1394  ! Concatenate offline directory and file names
1395  cs%snap_file = trim(cs%offlinedir)//trim(cs%snap_file)
1396  cs%mean_file = trim(cs%offlinedir)//trim(cs%mean_file)
1397  cs%sum_file = trim(cs%offlinedir)//trim(cs%sum_file)
1398  cs%surf_file = trim(cs%offlinedir)//trim(cs%surf_file)
1399 
1400  cs%num_vert_iter = cs%dt_offline/cs%dt_offline_vertical
1401 
1402  ! Map redistribute_method onto logicals in CS
1403  select case (redistribute_method)
1404  case ('barotropic')
1405  cs%redistribute_barotropic = .true.
1406  cs%redistribute_upwards = .false.
1407  case ('upwards')
1408  cs%redistribute_barotropic = .false.
1409  cs%redistribute_upwards = .true.
1410  case ('both')
1411  cs%redistribute_barotropic = .true.
1412  cs%redistribute_upwards = .true.
1413  case ('none')
1414  cs%redistribute_barotropic = .false.
1415  cs%redistribute_upwards = .false.
1416  end select
1417 
1418  ! Set the accumulated time to zero
1419  cs%accumulated_time = real_to_time(0.0)
1420  cs%vertical_time = cs%accumulated_time
1421  ! Set the starting read index for time-averaged and snapshotted fields
1422  cs%ridx_sum = cs%start_index
1423  if (cs%fields_are_offset) cs%ridx_snap = next_modulo_time(cs%start_index,cs%numtime)
1424  if (.not. cs%fields_are_offset) cs%ridx_snap = cs%start_index
1425 
1426  ! Copy members from other modules
1427  call extract_diabatic_member(diabatic_csp, opacity_csp=cs%opacity_CSp, optics_csp=cs%optics, &
1428  diabatic_aux_csp=cs%diabatic_aux_CSp, &
1429  evap_cfl_limit=cs%evap_CFL_limit, &
1430  minimum_forcing_depth=cs%minimum_forcing_depth)
1431 
1432  ! Grid pointer assignments
1433  cs%G => g
1434  cs%GV => gv
1435 
1436  ! Allocate arrays
1437  allocate(cs%uhtr(isdb:iedb,jsd:jed,nz)) ; cs%uhtr(:,:,:) = 0.0
1438  allocate(cs%vhtr(isd:ied,jsdb:jedb,nz)) ; cs%vhtr(:,:,:) = 0.0
1439  allocate(cs%eatr(isd:ied,jsd:jed,nz)) ; cs%eatr(:,:,:) = 0.0
1440  allocate(cs%ebtr(isd:ied,jsd:jed,nz)) ; cs%ebtr(:,:,:) = 0.0
1441  allocate(cs%h_end(isd:ied,jsd:jed,nz)) ; cs%h_end(:,:,:) = 0.0
1442  allocate(cs%netMassOut(g%isd:g%ied,g%jsd:g%jed)) ; cs%netMassOut(:,:) = 0.0
1443  allocate(cs%netMassIn(g%isd:g%ied,g%jsd:g%jed)) ; cs%netMassIn(:,:) = 0.0
1444  allocate(cs%Kd(isd:ied,jsd:jed,nz+1)) ; cs%Kd = 0.
1445  if (cs%read_mld) then
1446  allocate(cs%mld(g%isd:g%ied,g%jsd:g%jed)) ; cs%mld(:,:) = 0.0
1447  endif
1448 
1449  if (cs%read_all_ts_uvh) then
1450  call read_all_input(cs)
1451  endif
1452 
1453  ! Initialize ids for clocks used in offline routines
1454  cs%id_clock_read_fields = cpu_clock_id('(Offline read fields)',grain=clock_module)
1455  cs%id_clock_offline_diabatic = cpu_clock_id('(Offline diabatic)',grain=clock_module)
1456  cs%id_clock_offline_adv = cpu_clock_id('(Offline transport)',grain=clock_module)
1457  cs%id_clock_redistribute = cpu_clock_id('(Offline redistribute)',grain=clock_module)
1458 
1459  call calltree_leave("offline_transport_init")
1460 

◆ post_offline_convergence_diags()

subroutine, public mom_offline_main::post_offline_convergence_diags ( type(offline_transport_cs), intent(in)  CS,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  h_off,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  h_end,
real, dimension(szib_(cs%g),szj_(cs%g),szk_(cs%g)), intent(inout)  uhtr,
real, dimension(szi_(cs%g),szjb_(cs%g),szk_(cs%g)), intent(inout)  vhtr 
)

Posts diagnostics related to offline convergence diagnostics.

Parameters
[in]csOffline control structure
[in,out]h_offThicknesses at end of offline step
[in,out]h_endStored thicknesses
[in,out]uhtrRemaining zonal mass transport
[in,out]vhtrRemaining meridional mass transport

Definition at line 1172 of file MOM_offline_main.F90.

1172  type(offline_transport_CS), intent(in ) :: CS !< Offline control structure
1173  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_off !< Thicknesses at end of offline step
1174  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: h_end !< Stored thicknesses
1175  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(inout) :: uhtr !< Remaining zonal mass transport
1176  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(inout) :: vhtr !< Remaining meridional mass transport
1177 
1178  real, dimension(SZI_(CS%G),SZJ_(CS%G)) :: eta_diff
1179  integer :: i, j, k
1180 
1181  if (cs%id_eta_diff_end>0) then
1182  ! Calculate difference in column thickness
1183  eta_diff = 0.
1184  do k=1,cs%GV%ke ; do j=cs%G%jsc,cs%G%jec ; do i=cs%G%isc,cs%G%iec
1185  eta_diff(i,j) = eta_diff(i,j) + h_off(i,j,k)
1186  enddo ; enddo ; enddo
1187  do k=1,cs%GV%ke ; do j=cs%G%jsc,cs%G%jec ; do i=cs%G%isc,cs%G%iec
1188  eta_diff(i,j) = eta_diff(i,j) - h_end(i,j,k)
1189  enddo ; enddo ; enddo
1190 
1191  call post_data(cs%id_eta_diff_end, eta_diff, cs%diag)
1192  endif
1193 
1194  if (cs%id_hdiff>0) call post_data(cs%id_hdiff, h_off-h_end, cs%diag)
1195  if (cs%id_hr>0) call post_data(cs%id_hr, h_off, cs%diag)
1196  if (cs%id_uhr_end>0) call post_data(cs%id_uhr_end, uhtr, cs%diag)
1197  if (cs%id_vhr_end>0) call post_data(cs%id_vhr_end, vhtr, cs%diag)
1198 

◆ read_all_input()

subroutine mom_offline_main::read_all_input ( type(offline_transport_cs), intent(inout)  CS)
private

Coordinates the allocation and reading in all time levels of uh, vh, hend, temp, and salt from files. Used when read_all_ts_uvh.

Parameters
[in,out]csControl structure for offline module

Definition at line 1466 of file MOM_offline_main.F90.

1466  type(offline_transport_CS), intent(inout) :: CS !< Control structure for offline module
1467 
1468  integer :: is, ie, js, je, isd, ied, jsd, jed, nz, t, ntime
1469  integer :: IsdB, IedB, JsdB, JedB
1470 
1471  nz = cs%GV%ke ; ntime = cs%numtime
1472  isd = cs%G%isd ; ied = cs%G%ied ; jsd = cs%G%jsd ; jed = cs%G%jed
1473  isdb = cs%G%IsdB ; iedb = cs%G%IedB ; jsdb = cs%G%JsdB ; jedb = cs%G%JedB
1474 
1475  ! Extra safety check that we're not going to overallocate any arrays
1476  if (cs%read_all_ts_uvh) then
1477  if (allocated(cs%uhtr_all)) call mom_error(fatal, "uhtr_all is already allocated")
1478  if (allocated(cs%vhtr_all)) call mom_error(fatal, "vhtr_all is already allocated")
1479  if (allocated(cs%hend_all)) call mom_error(fatal, "hend_all is already allocated")
1480  if (allocated(cs%temp_all)) call mom_error(fatal, "temp_all is already allocated")
1481  if (allocated(cs%salt_all)) call mom_error(fatal, "salt_all is already allocated")
1482 
1483  allocate(cs%uhtr_all(isdb:iedb,jsd:jed,nz,ntime)) ; cs%uhtr_all(:,:,:,:) = 0.0
1484  allocate(cs%vhtr_all(isd:ied,jsdb:jedb,nz,ntime)) ; cs%vhtr_all(:,:,:,:) = 0.0
1485  allocate(cs%hend_all(isd:ied,jsd:jed,nz,ntime)) ; cs%hend_all(:,:,:,:) = 0.0
1486  allocate(cs%temp_all(isd:ied,jsd:jed,nz,1:ntime)) ; cs%temp_all(:,:,:,:) = 0.0
1487  allocate(cs%salt_all(isd:ied,jsd:jed,nz,1:ntime)) ; cs%salt_all(:,:,:,:) = 0.0
1488 
1489  call mom_mesg("Reading in uhtr, vhtr, h_start, h_end, temp, salt")
1490  do t = 1,ntime
1491  call mom_read_vector(cs%snap_file, 'uhtr_sum', 'vhtr_sum', cs%uhtr_all(:,:,1:cs%nk_input,t), &
1492  cs%vhtr_all(:,:,1:cs%nk_input,t), cs%G%Domain, timelevel=t)
1493  call mom_read_data(cs%snap_file,'h_end', cs%hend_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1494  timelevel=t, position=center)
1495  call mom_read_data(cs%mean_file,'temp', cs%temp_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1496  timelevel=t, position=center)
1497  call mom_read_data(cs%mean_file,'salt', cs%salt_all(:,:,1:cs%nk_input,t), cs%G%Domain, &
1498  timelevel=t, position=center)
1499  enddo
1500  endif
1501 

◆ register_diags_offline_transport()

subroutine, public mom_offline_main::register_diags_offline_transport ( type(time_type), intent(in)  Time,
type(diag_ctrl), intent(in)  diag,
type(offline_transport_cs), pointer  CS 
)

Initialize additional diagnostics required for offline tracer transport.

Parameters
csControl structure for offline module
[in]timecurrent model time
[in]diagStructure that regulates diagnostic output

Definition at line 1115 of file MOM_offline_main.F90.

1115 
1116  type(offline_transport_CS), pointer :: CS !< Control structure for offline module
1117  type(time_type), intent(in) :: Time !< current model time
1118  type(diag_ctrl), intent(in) :: diag !< Structure that regulates diagnostic output
1119 
1120  ! U-cell fields
1121  cs%id_uhr = register_diag_field('ocean_model', 'uhr', diag%axesCuL, time, &
1122  'Zonal thickness fluxes remaining at end of advection', 'kg')
1123  cs%id_uhr_redist = register_diag_field('ocean_model', 'uhr_redist', diag%axesCuL, time, &
1124  'Zonal thickness fluxes to be redistributed vertically', 'kg')
1125  cs%id_uhr_end = register_diag_field('ocean_model', 'uhr_end', diag%axesCuL, time, &
1126  'Zonal thickness fluxes at end of offline step', 'kg')
1127 
1128  ! V-cell fields
1129  cs%id_vhr = register_diag_field('ocean_model', 'vhr', diag%axesCvL, time, &
1130  'Meridional thickness fluxes remaining at end of advection', 'kg')
1131  cs%id_vhr_redist = register_diag_field('ocean_model', 'vhr_redist', diag%axesCvL, time, &
1132  'Meridional thickness to be redistributed vertically', 'kg')
1133  cs%id_vhr_end = register_diag_field('ocean_model', 'vhr_end', diag%axesCvL, time, &
1134  'Meridional thickness at end of offline step', 'kg')
1135 
1136  ! T-cell fields
1137  cs%id_hdiff = register_diag_field('ocean_model', 'hdiff', diag%axesTL, time, &
1138  'Difference between the stored and calculated layer thickness', 'm')
1139  cs%id_hr = register_diag_field('ocean_model', 'hr', diag%axesTL, time, &
1140  'Layer thickness at end of offline step', 'm')
1141  cs%id_ear = register_diag_field('ocean_model', 'ear', diag%axesTL, time, &
1142  'Remaining thickness entrained from above', 'm')
1143  cs%id_ebr = register_diag_field('ocean_model', 'ebr', diag%axesTL, time, &
1144  'Remaining thickness entrained from below', 'm')
1145  cs%id_eta_pre_distribute = register_diag_field('ocean_model','eta_pre_distribute', &
1146  diag%axesT1, time, 'Total water column height before residual transport redistribution','m')
1147  cs%id_eta_post_distribute = register_diag_field('ocean_model','eta_post_distribute', &
1148  diag%axesT1, time, 'Total water column height after residual transport redistribution','m')
1149  cs%id_eta_diff_end = register_diag_field('ocean_model','eta_diff_end', diag%axesT1, time, &
1150  'Difference in total water column height from online and offline ' // &
1151  'at the end of the offline timestep','m')
1152  cs%id_h_redist = register_diag_field('ocean_model','h_redist', diag%axesTL, time, &
1153  'Layer thicknesses before redistribution of mass fluxes','m')
1154 
1155  ! Regridded/remapped input fields
1156  cs%id_uhtr_regrid = register_diag_field('ocean_model', 'uhtr_regrid', diag%axesCuL, time, &
1157  'Zonal mass transport regridded/remapped onto offline grid','kg')
1158  cs%id_vhtr_regrid = register_diag_field('ocean_model', 'vhtr_regrid', diag%axesCvL, time, &
1159  'Meridional mass transport regridded/remapped onto offline grid','kg')
1160  cs%id_temp_regrid = register_diag_field('ocean_model', 'temp_regrid', diag%axesTL, time, &
1161  'Temperature regridded/remapped onto offline grid','C')
1162  cs%id_salt_regrid = register_diag_field('ocean_model', 'salt_regrid', diag%axesTL, time, &
1163  'Salinity regridded/remapped onto offline grid','g kg-1')
1164  cs%id_h_regrid = register_diag_field('ocean_model', 'h_regrid', diag%axesTL, time, &
1165  'Layer thicknesses regridded/remapped onto offline grid','m')
1166 
1167 

◆ remaining_transport_sum()

real function mom_offline_main::remaining_transport_sum ( type(offline_transport_cs), pointer  CS,
real, dimension(szib_(cs%g),szj_(cs%g),szk_(cs%g)), intent(in)  uhtr,
real, dimension(szi_(cs%g),szjb_(cs%g),szk_(cs%g)), intent(in)  vhtr 
)
private

Sums any non-negligible remaining transport to check for advection convergence.

Parameters
cscontrol structure for offline module
[in]uhtrZonal mass transport
[in]vhtrMeridional mass transport

Definition at line 619 of file MOM_offline_main.F90.

619  type(offline_transport_CS), pointer :: CS !< control structure for offline module
620  real, dimension(SZIB_(CS%G),SZJ_(CS%G),SZK_(CS%G)), intent(in ) :: uhtr !< Zonal mass transport
621  real, dimension(SZI_(CS%G),SZJB_(CS%G),SZK_(CS%G)), intent(in ) :: vhtr !< Meridional mass transport
622 
623  ! Local variables
624  integer :: i, j, k
625  integer :: is, ie, js, je, nz
626  real :: h_min !< A layer thickness below roundoff from GV type
627  real :: uh_neglect !< A small value of zonal transport that effectively is below roundoff error
628  real :: vh_neglect !< A small value of meridional transport that effectively is below roundoff error
629 
630  nz = cs%GV%ke
631  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec
632 
633  h_min = cs%GV%H_subroundoff
634 
635  remaining_transport_sum = 0.
636  do k=1,nz; do j=js,je ; do i=is,ie
637  uh_neglect = h_min*cs%G%US%L_to_m**2*min(cs%G%areaT(i,j),cs%G%areaT(i+1,j))
638  vh_neglect = h_min*cs%G%US%L_to_m**2*min(cs%G%areaT(i,j),cs%G%areaT(i,j+1))
639  if (abs(uhtr(i,j,k))>uh_neglect) then
640  remaining_transport_sum = remaining_transport_sum + abs(uhtr(i,j,k))
641  endif
642  if (abs(vhtr(i,j,k))>vh_neglect) then
643  remaining_transport_sum = remaining_transport_sum + abs(vhtr(i,j,k))
644  endif
645  enddo ; enddo ; enddo
646  call sum_across_pes(remaining_transport_sum)
647 

◆ update_offline_fields()

subroutine, public mom_offline_main::update_offline_fields ( type(offline_transport_cs), pointer  CS,
real, dimension(szi_(cs%g),szj_(cs%g),szk_(cs%g))  h,
type(forcing), intent(inout)  fluxes,
logical, intent(in)  do_ale 
)

Update fields used in this round of offline transport. First fields are updated from files or from arrays read during initialization. Then if in an ALE-dependent coordinate, regrid/remap fields.

Parameters
csControl structure for offline module
hThe regridded layer thicknesses
[in,out]fluxesPointers to forcing fields
[in]do_aleTrue if using ALE

Definition at line 1019 of file MOM_offline_main.F90.

1019  type(offline_transport_CS), pointer :: CS !< Control structure for offline module
1020  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: h !< The regridded layer thicknesses
1021  type(forcing), intent(inout) :: fluxes !< Pointers to forcing fields
1022  logical, intent(in ) :: do_ale !< True if using ALE
1023  ! Local variables
1024  integer :: i, j, k, is, ie, js, je, nz
1025  real, dimension(SZI_(CS%G),SZJ_(CS%G),SZK_(CS%G)) :: h_start
1026  is = cs%G%isc ; ie = cs%G%iec ; js = cs%G%jsc ; je = cs%G%jec ; nz = cs%GV%ke
1027 
1028  call cpu_clock_begin(cs%id_clock_read_fields)
1029  call calltree_enter("update_offline_fields, MOM_offline_main.F90")
1030 
1031  ! Store a copy of the layer thicknesses before ALE regrid/remap
1032  h_start(:,:,:) = h(:,:,:)
1033 
1034  ! Most fields will be read in from files
1035  call update_offline_from_files( cs%G, cs%GV, cs%nk_input, cs%mean_file, cs%sum_file, cs%snap_file, cs%surf_file, &
1036  cs%h_end, cs%uhtr, cs%vhtr, cs%tv%T, cs%tv%S, cs%mld, cs%Kd, fluxes, &
1037  cs%ridx_sum, cs%ridx_snap, cs%read_mld, cs%read_sw, .not. cs%read_all_ts_uvh, do_ale)
1038  ! If uh, vh, h_end, temp, salt were read in at the beginning, fields are copied from those arrays
1039  if (cs%read_all_ts_uvh) then
1040  call update_offline_from_arrays(cs%G, cs%GV, cs%nk_input, cs%ridx_sum, cs%mean_file, cs%sum_file, cs%snap_file, &
1041  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)
1042  endif
1043  if (cs%debug) then
1044  call uvchksum("[uv]h after update offline from files and arrays", cs%uhtr, cs%vhtr, cs%G%HI)
1045  endif
1046 
1047  ! If using an ALE-dependent vertical coordinate, fields will need to be remapped
1048  if (do_ale) then
1049  ! These halo passes are necessary because u, v fields will need information 1 step into the halo
1050  call pass_var(h, cs%G%Domain)
1051  call pass_var(cs%tv%T, cs%G%Domain)
1052  call pass_var(cs%tv%S, cs%G%Domain)
1053  call ale_offline_inputs(cs%ALE_CSp, cs%G, cs%GV, h, cs%tv, cs%tracer_Reg, cs%uhtr, cs%vhtr, cs%Kd, &
1054  cs%debug, cs%OBC)
1055  if (cs%id_temp_regrid>0) call post_data(cs%id_temp_regrid, cs%tv%T, cs%diag)
1056  if (cs%id_salt_regrid>0) call post_data(cs%id_salt_regrid, cs%tv%S, cs%diag)
1057  if (cs%id_uhtr_regrid>0) call post_data(cs%id_uhtr_regrid, cs%uhtr, cs%diag)
1058  if (cs%id_vhtr_regrid>0) call post_data(cs%id_vhtr_regrid, cs%vhtr, cs%diag)
1059  if (cs%id_h_regrid>0) call post_data(cs%id_h_regrid, h, cs%diag)
1060  if (cs%debug) then
1061  call uvchksum("[uv]h after ALE regridding/remapping of inputs", cs%uhtr, cs%vhtr, cs%G%HI)
1062  call hchksum(h_start,"h_start after update offline from files and arrays", cs%G%HI)
1063  endif
1064  endif
1065 
1066  ! Update halos for some
1067  call pass_var(cs%h_end, cs%G%Domain)
1068  call pass_var(cs%tv%T, cs%G%Domain)
1069  call pass_var(cs%tv%S, cs%G%Domain)
1070 
1071  ! Update the read indices
1072  cs%ridx_snap = next_modulo_time(cs%ridx_snap,cs%numtime)
1073  cs%ridx_sum = next_modulo_time(cs%ridx_sum,cs%numtime)
1074 
1075  ! Apply masks/factors at T, U, and V points
1076  do k=1,nz ; do j=js,je ; do i=is,ie
1077  if (cs%G%mask2dT(i,j)<1.0) then
1078  cs%h_end(i,j,k) = cs%GV%Angstrom_H
1079  endif
1080  enddo ; enddo ; enddo
1081 
1082  do k=1,nz+1 ; do j=js,je ; do i=is,ie
1083  cs%Kd(i,j,k) = max(0.0, cs%Kd(i,j,k))
1084  if (cs%Kd_max>0.) then
1085  cs%Kd(i,j,k) = min(cs%Kd_max, cs%Kd(i,j,k))
1086  endif
1087  enddo ; enddo ; enddo
1088 
1089  do k=1,nz ; do j=js-1,je ; do i=is,ie
1090  if (cs%G%mask2dCv(i,j)<1.0) then
1091  cs%vhtr(i,j,k) = 0.0
1092  endif
1093  enddo ; enddo ; enddo
1094 
1095  do k=1,nz ; do j=js,je ; do i=is-1,ie
1096  if (cs%G%mask2dCu(i,j)<1.0) then
1097  cs%uhtr(i,j,k) = 0.0
1098  endif
1099  enddo ; enddo ; enddo
1100 
1101  if (cs%debug) then
1102  call uvchksum("[uv]htr_sub after update_offline_fields", cs%uhtr, cs%vhtr, cs%G%HI)
1103  call hchksum(cs%h_end, "h_end after update_offline_fields", cs%G%HI)
1104  call hchksum(cs%tv%T, "Temp after update_offline_fields", cs%G%HI)
1105  call hchksum(cs%tv%S, "Salt after update_offline_fields", cs%G%HI)
1106  endif
1107 
1108  call calltree_leave("update_offline_fields")
1109  call cpu_clock_end(cs%id_clock_read_fields)
1110