MOM6
MOM_offline_aux.F90
1 !> Contains routines related to offline transport of tracers. These routines are likely to be called from
2 !> the MOM_offline_main module
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use data_override_mod, only : data_override_init, data_override
8 use mom_time_manager, only : time_type, operator(-)
9 use mom_debugging, only : check_column_integrals
10 use mom_domains, only : pass_var, pass_vector, to_all
11 use mom_diag_vkernels, only : reintegrate_column
12 use mom_error_handler, only : calltree_enter, calltree_leave, mom_error, fatal, warning, is_root_pe
13 use mom_grid, only : ocean_grid_type
14 use mom_io, only : mom_read_data, mom_read_vector, center
17 use astronomy_mod, only : orbital_time, diurnal_solar, daily_mean_solar
18 use mom_variables, only : vertvisc_type
19 use mom_forcing_type, only : forcing
20 use mom_opacity, only : optics_type
21 use mom_diag_mediator, only : post_data
22 use mom_forcing_type, only : forcing
23 
24 implicit none ; private
25 
26 public update_offline_from_files
27 public update_offline_from_arrays
28 public update_h_horizontal_flux
29 public update_h_vertical_flux
30 public limit_mass_flux_3d
31 public distribute_residual_uh_barotropic
32 public distribute_residual_vh_barotropic
33 public distribute_residual_uh_upwards
34 public distribute_residual_vh_upwards
35 public next_modulo_time
36 public offline_add_diurnal_sw
37 
38 #include "MOM_memory.h"
39 #include "version_variable.h"
40 
41 contains
42 
43 !> This updates thickness based on the convergence of horizontal mass fluxes
44 !! NOTE: Only used in non-ALE mode
45 subroutine update_h_horizontal_flux(G, GV, uhtr, vhtr, h_pre, h_new)
46  type(ocean_grid_type), pointer :: g !< ocean grid structure
47  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
48  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
49  intent(in) :: uhtr !< Accumulated mass flux through zonal face [kg]
50  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
51  intent(in) :: vhtr !< Accumulated mass flux through meridional face [kg]
52  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
53  intent(in) :: h_pre !< Previous layer thicknesses [kg m-2].
54  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
55  intent(inout) :: h_new !< Updated layer thicknesses [kg m-2].
56 
57  ! Local variables
58  integer :: i, j, k, m, is, ie, js, je, nz
59  ! Set index-related variables for fields on T-grid
60  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
61 
62  do k = 1, nz
63  do i=is-1,ie+1 ; do j=js-1,je+1
64 
65  h_new(i,j,k) = max(0.0, g%US%L_to_m**2*g%areaT(i,j)*h_pre(i,j,k) + &
66  ((uhtr(i-1,j,k) - uhtr(i,j,k)) + (vhtr(i,j-1,k) - vhtr(i,j,k))))
67 
68  ! In the case that the layer is now dramatically thinner than it was previously,
69  ! add a bit of mass to avoid truncation errors. This will lead to
70  ! non-conservation of tracers
71  h_new(i,j,k) = h_new(i,j,k) + &
72  max(gv%Angstrom_H, 1.0e-13*h_new(i,j,k) - g%US%L_to_m**2*g%areaT(i,j)*h_pre(i,j,k))
73 
74  ! Convert back to thickness
75  h_new(i,j,k) = h_new(i,j,k) / (g%US%L_to_m**2*g%areaT(i,j))
76 
77  enddo ; enddo
78  enddo
79 end subroutine update_h_horizontal_flux
80 
81 !> Updates layer thicknesses due to vertical mass transports
82 !! NOTE: Only used in non-ALE configuration
83 subroutine update_h_vertical_flux(G, GV, ea, eb, h_pre, h_new)
84  type(ocean_grid_type), pointer :: g !< ocean grid structure
85  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
86  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
87  intent(in) :: ea !< Mass of fluid entrained from the layer
88  !! above within this timestep [kg m-2]
89  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
90  intent(in) :: eb !< Mass of fluid entrained from the layer
91  !! below within this timestep [kg m-2]
92  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
93  intent(in) :: h_pre !< Layer thicknesses at the end of the previous
94  !! step [kg m-2].
95  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
96  intent(inout) :: h_new !< Updated layer thicknesses [kg m-2].
97 
98  ! Local variables
99  integer :: i, j, k, m, is, ie, js, je, nz
100  ! Set index-related variables for fields on T-grid
101  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
102 
103  ! Update h_new with convergence of vertical mass transports
104  do j=js-1,je+1
105  do i=is-1,ie+1
106 
107  ! Top layer
108  h_new(i,j,1) = max(0.0, h_pre(i,j,1) + (eb(i,j,1) - ea(i,j,2) + ea(i,j,1) ))
109  h_new(i,j,1) = h_new(i,j,1) + &
110  max(0.0, 1.0e-13*h_new(i,j,1) - h_pre(i,j,1))
111 
112  ! Bottom layer
113 ! h_new(i,j,nz) = h_pre(i,j,nz) + (ea(i,j,nz) - eb(i,j,nz-1)+eb(i,j,nz))
114  h_new(i,j,nz) = max(0.0, h_pre(i,j,nz) + (ea(i,j,nz) - eb(i,j,nz-1)+eb(i,j,nz)))
115  h_new(i,j,nz) = h_new(i,j,nz) + &
116  max(0.0, 1.0e-13*h_new(i,j,nz) - h_pre(i,j,nz))
117 
118  enddo
119 
120  ! Interior layers
121  do k=2,nz-1 ; do i=is-1,ie+1
122 
123  h_new(i,j,k) = max(0.0, h_pre(i,j,k) + ((ea(i,j,k) - eb(i,j,k-1)) + &
124  (eb(i,j,k) - ea(i,j,k+1))))
125  h_new(i,j,k) = h_new(i,j,k) + &
126  max(0.0, 1.0e-13*h_new(i,j,k) - h_pre(i,j,k))
127 
128  enddo ; enddo
129 
130  enddo
131 
132 end subroutine update_h_vertical_flux
133 
134 !> This routine limits the mass fluxes so that the a layer cannot be completely depleted.
135 !! NOTE: Only used in non-ALE mode
136 subroutine limit_mass_flux_3d(G, GV, uh, vh, ea, eb, h_pre)
137  type(ocean_grid_type), pointer :: g !< ocean grid structure
138  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
139  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
140  intent(inout) :: uh !< Mass flux through zonal face [kg]
141  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
142  intent(inout) :: vh !< Mass flux through meridional face [kg]
143  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
144  intent(inout) :: ea !< Mass of fluid entrained from the layer
145  !! above within this timestep [kg m-2]
146  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
147  intent(inout) :: eb !< Mass of fluid entrained from the layer
148  !! below within this timestep [kg m-2]
149  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
150  intent(in) :: h_pre !< Layer thicknesses at the end of the previous
151  !! step [kg m-2].
152 
153  ! Local variables
154  integer :: i, j, k, m, is, ie, js, je, nz
155  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: top_flux, bottom_flux
156  real :: pos_flux, hvol, h_neglect, scale_factor, max_off_cfl
157 
158  max_off_cfl =0.5
159 
160  ! In this subroutine, fluxes out of the box are scaled away if they deplete
161  ! the layer, note that we define the positive direction as flux out of the box.
162  ! Hence, uh(I-1) is multipled by negative one, but uh(I) is not
163 
164  ! Set index-related variables for fields on T-grid
165  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
166 
167  ! Calculate top and bottom fluxes from ea and eb. Note the explicit negative signs
168  ! to enforce the positive out convention
169  k = 1
170  do j=js-1,je+1 ; do i=is-1,ie+1
171  top_flux(i,j,k) = -ea(i,j,k)
172  bottom_flux(i,j,k) = -(eb(i,j,k)-ea(i,j,k+1))
173  enddo ; enddo
174 
175  do k=2, nz-1 ; do j=js-1,je+1 ; do i=is-1,ie+1
176  top_flux(i,j,k) = -(ea(i,j,k)-eb(i,j,k-1))
177  bottom_flux(i,j,k) = -(eb(i,j,k)-ea(i,j,k+1))
178  enddo ; enddo ; enddo
179 
180  k=nz
181  do j=js-1,je+1 ; do i=is-1,ie+1
182  top_flux(i,j,k) = -(ea(i,j,k)-eb(i,j,k-1))
183  bottom_flux(i,j,k) = -eb(i,j,k)
184  enddo ; enddo
185 
186 
187  ! Calculate sum of positive fluxes (negatives applied to enforce convention)
188  ! in a given cell and scale it back if it would deplete a layer
189  do k = 1, nz ; do j=js-1,je+1 ; do i=is-1,ie+1
190 
191  hvol = h_pre(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)
192  pos_flux = max(0.0,-uh(i-1,j,k)) + max(0.0, -vh(i,j-1,k)) + &
193  max(0.0, uh(i,j,k)) + max(0.0, vh(i,j,k)) + &
194  max(0.0, top_flux(i,j,k)*g%US%L_to_m**2*g%areaT(i,j)) + max(0.0, bottom_flux(i,j,k)*g%US%L_to_m**2*g%areaT(i,j))
195 
196  if (pos_flux>hvol .and. pos_flux>0.0) then
197  scale_factor = ( hvol )/pos_flux*max_off_cfl
198  else ! Don't scale
199  scale_factor = 1.0
200  endif
201 
202  ! Scale horizontal fluxes
203  if (-uh(i-1,j,k)>0) uh(i-1,j,k) = uh(i-1,j,k)*scale_factor
204  if (uh(i,j,k)>0) uh(i,j,k) = uh(i,j,k)*scale_factor
205  if (-vh(i,j-1,k)>0) vh(i,j-1,k) = vh(i,j-1,k)*scale_factor
206  if (vh(i,j,k)>0) vh(i,j,k) = vh(i,j,k)*scale_factor
207 
208  if (k>1 .and. k<nz) then
209  ! Scale interior layers
210  if (top_flux(i,j,k)>0.0) then
211  ea(i,j,k) = ea(i,j,k)*scale_factor
212  eb(i,j,k-1) = eb(i,j,k-1)*scale_factor
213  endif
214  if (bottom_flux(i,j,k)>0.0) then
215  eb(i,j,k) = eb(i,j,k)*scale_factor
216  ea(i,j,k+1) = ea(i,j,k+1)*scale_factor
217  endif
218  ! Scale top layer
219  elseif (k==1) then
220  if (top_flux(i,j,k)>0.0) ea(i,j,k) = ea(i,j,k)*scale_factor
221  if (bottom_flux(i,j,k)>0.0) then
222  eb(i,j,k) = eb(i,j,k)*scale_factor
223  ea(i,j,k+1) = ea(i,j,k+1)*scale_factor
224  endif
225  ! Scale bottom layer
226  elseif (k==nz) then
227  if (top_flux(i,j,k)>0.0) then
228  ea(i,j,k) = ea(i,j,k)*scale_factor
229  eb(i,j,k-1) = eb(i,j,k-1)*scale_factor
230  endif
231  if (bottom_flux(i,j,k)>0.0) eb(i,j,k)=eb(i,j,k)*scale_factor
232  endif
233  enddo ; enddo ; enddo
234 
235 end subroutine limit_mass_flux_3d
236 
237 !> In the case where offline advection has failed to converge, redistribute the u-flux
238 !! into remainder of the water column as a barotropic equivalent
239 subroutine distribute_residual_uh_barotropic(G, GV, hvol, uh)
240  type(ocean_grid_type), pointer :: g !< ocean grid structure
241  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
242  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
243  intent(in ) :: hvol !< Mass of water in the cells at the end
244  !! of the previous timestep [kg]
245  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
246  intent(inout) :: uh !< Zonal mass transport within a timestep [kg]
247 
248  real, dimension(SZIB_(G),SZK_(G)) :: uh2d
249  real, dimension(SZIB_(G)) :: uh2d_sum
250  real, dimension(SZI_(G),SZK_(G)) :: h2d
251  real, dimension(SZI_(G)) :: h2d_sum
252 
253  integer :: i, j, k, m, is, ie, js, je, nz
254  real :: uh_neglect
255 
256  ! Set index-related variables for fields on T-grid
257  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
258 
259  do j=js,je
260  uh2d_sum(:) = 0.0
261  ! Copy over uh to a working array and sum up the remaining fluxes in a column
262  do k=1,nz ; do i=is-1,ie
263  uh2d(i,k) = uh(i,j,k)
264  uh2d_sum(i) = uh2d_sum(i) + uh2d(i,k)
265  enddo ; enddo
266 
267  ! Copy over h to a working array and calculate total column volume
268  h2d_sum(:) = 0.0
269  do k=1,nz ; do i=is-1,ie+1
270  h2d(i,k) = hvol(i,j,k)
271  if (hvol(i,j,k)>0.) then
272  h2d_sum(i) = h2d_sum(i) + h2d(i,k)
273  else
274  h2d(i,k) = gv%H_subroundoff
275  endif
276  enddo ; enddo
277 
278  ! Distribute flux. Note min/max is intended to make sure that the mass transport
279  ! does not deplete a cell
280  do i=is-1,ie
281  if ( uh2d_sum(i)>0.0 ) then
282  do k=1,nz
283  uh2d(i,k) = uh2d_sum(i)*(h2d(i,k)/h2d_sum(i))
284  enddo
285  elseif (uh2d_sum(i)<0.0) then
286  do k=1,nz
287  uh2d(i,k) = uh2d_sum(i)*(h2d(i+1,k)/h2d_sum(i+1))
288  enddo
289  else
290  do k=1,nz
291  uh2d(i,k) = 0.0
292  enddo
293  endif
294  ! Calculate and check that column integrated transports match the original to
295  ! within the tolerance limit
296  uh_neglect = gv%Angstrom_H*g%US%L_to_m**2 * min(g%areaT(i,j), g%areaT(i+1,j))
297  if ( abs(sum(uh2d(i,:))-uh2d_sum(i)) > uh_neglect) &
298  call mom_error(warning,"Column integral of uh does not match after "//&
299  "barotropic redistribution")
300  enddo
301 
302  do k=1,nz ; do i=is-1,ie
303  uh(i,j,k) = uh2d(i,k)
304  enddo ; enddo
305  enddo
306 
307 end subroutine distribute_residual_uh_barotropic
308 
309 !> Redistribute the v-flux as a barotropic equivalent
310 subroutine distribute_residual_vh_barotropic(G, GV, hvol, vh)
311  type(ocean_grid_type), pointer :: g !< ocean grid structure
312  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
313  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
314  intent(in ) :: hvol !< Mass of water in the cells at the end
315  !! of the previous timestep [kg]
316  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
317  intent(inout) :: vh !< Meridional mass transport within a timestep [kg]
318 
319  real, dimension(SZJB_(G),SZK_(G)) :: vh2d
320  real, dimension(SZJB_(G)) :: vh2d_sum
321  real, dimension(SZJ_(G),SZK_(G)) :: h2d
322  real, dimension(SZJ_(G)) :: h2d_sum
323 
324  integer :: i, j, k, m, is, ie, js, je, nz
325  real :: vh_neglect
326 
327  ! Set index-related variables for fields on T-grid
328  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
329 
330  do i=is,ie
331  vh2d_sum(:) = 0.0
332  ! Copy over uh to a working array and sum up the remaining fluxes in a column
333  do k=1,nz ; do j=js-1,je
334  vh2d(j,k) = vh(i,j,k)
335  vh2d_sum(j) = vh2d_sum(j) + vh2d(j,k)
336  enddo ; enddo
337 
338  ! Copy over h to a working array and calculate column volume
339  h2d_sum(:) = 0.0
340  do k=1,nz ; do j=js-1,je+1
341  h2d(j,k) = hvol(i,j,k)
342  if (hvol(i,j,k)>0.) then
343  h2d_sum(j) = h2d_sum(j) + h2d(j,k)
344  else
345  h2d(j,k) = gv%H_subroundoff
346  endif
347  enddo ; enddo
348 
349  ! Distribute flux evenly throughout a column
350  do j=js-1,je
351  if ( vh2d_sum(j)>0.0 ) then
352  do k=1,nz
353  vh2d(j,k) = vh2d_sum(j)*(h2d(j,k)/h2d_sum(j))
354  enddo
355  elseif (vh2d_sum(j)<0.0) then
356  do k=1,nz
357  vh2d(j,k) = vh2d_sum(j)*(h2d(j+1,k)/h2d_sum(j+1))
358  enddo
359  else
360  do k=1,nz
361  vh2d(j,k) = 0.0
362  enddo
363  endif
364  ! Calculate and check that column integrated transports match the original to
365  ! within the tolerance limit
366  vh_neglect = gv%Angstrom_H*g%US%L_to_m**2 * min(g%areaT(i,j), g%areaT(i,j+1))
367  if ( abs(sum(vh2d(j,:))-vh2d_sum(j)) > vh_neglect) then
368  call mom_error(warning,"Column integral of vh does not match after "//&
369  "barotropic redistribution")
370  endif
371 
372  enddo
373 
374  do k=1,nz ; do j=js-1,je
375  vh(i,j,k) = vh2d(j,k)
376  enddo ; enddo
377  enddo
378 
379 end subroutine distribute_residual_vh_barotropic
380 
381 !> In the case where offline advection has failed to converge, redistribute the u-flux
382 !! into layers above
383 subroutine distribute_residual_uh_upwards(G, GV, hvol, uh)
384  type(ocean_grid_type), pointer :: g !< ocean grid structure
385  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
386  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
387  intent(in ) :: hvol !< Mass of water in the cells at the end
388  !! of the previous timestep [kg]
389  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
390  intent(inout) :: uh !< Zonal mass transport within a timestep [kg]
391 
392  real, dimension(SZIB_(G),SZK_(G)) :: uh2d
393  real, dimension(SZI_(G),SZK_(G)) :: h2d
394 
395  real :: uh_neglect, uh_remain, uh_add, uh_sum, uh_col, uh_max
396  real :: hup, hdown, hlos, min_h
397  integer :: i, j, k, m, is, ie, js, je, nz, k_rev
398 
399  ! Set index-related variables for fields on T-grid
400  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
401 
402  min_h = gv%Angstrom_H*0.1
403 
404  do j=js,je
405  ! Copy over uh and cell volume to working arrays
406  do k=1,nz ; do i=is-2,ie+1
407  uh2d(i,k) = uh(i,j,k)
408  enddo ; enddo
409  do k=1,nz ; do i=is-1,ie+1
410  ! Subtract just a little bit of thickness to avoid roundoff errors
411  h2d(i,k) = hvol(i,j,k)-min_h*g%US%L_to_m**2*g%areaT(i,j)
412  enddo ; enddo
413 
414  do i=is-1,ie
415  uh_col = sum(uh2d(i,:)) ! Store original column-integrated transport
416  do k=1,nz
417  uh_remain = uh2d(i,k)
418  uh2d(i,k) = 0.0
419  if (abs(uh_remain)>0.0) then
420  do k_rev = k,1,-1
421  uh_sum = uh_remain + uh2d(i,k_rev)
422  if (uh_sum<0.0) then ! Transport to the left
423  hup = h2d(i+1,k_rev)
424  hlos = max(0.0,uh2d(i+1,k_rev))
425  if ((((hup - hlos) + uh_sum) < 0.0) .and. &
426  ((0.5*hup + uh_sum) < 0.0)) then
427  uh2d(i,k_rev) = min(-0.5*hup,-hup+hlos,0.0)
428  uh_remain = uh_sum - uh2d(i,k_rev)
429  else
430  uh2d(i,k_rev) = uh_sum
431  uh_remain = 0.0
432  exit
433  endif
434  else ! Transport to the right
435  hup = h2d(i,k_rev)
436  hlos = max(0.0,-uh2d(i-1,k_rev))
437  if ((((hup - hlos) - uh_sum) < 0.0) .and. &
438  ((0.5*hup - uh_sum) < 0.0)) then
439  uh2d(i,k_rev) = max(0.5*hup,hup-hlos,0.0)
440  uh_remain = uh_sum - uh2d(i,k_rev)
441  else
442  uh2d(i,k_rev) = uh_sum
443  uh_remain = 0.0
444  exit
445  endif
446  endif
447  enddo ! k_rev
448  endif
449 
450  if (abs(uh_remain)>0.0) then
451  if (k<nz) then
452  uh2d(i,k+1) = uh2d(i,k+1) + uh_remain
453  else
454  uh2d(i,k) = uh2d(i,k) + uh_remain
455  call mom_error(warning,"Water column cannot accommodate UH redistribution. Tracer may not be conserved")
456  endif
457  endif
458  enddo ! k-loop
459 
460  ! Calculate and check that column integrated transports match the original to
461  ! within the tolerance limit
462  uh_neglect = gv%Angstrom_H*g%US%L_to_m**2 * min(g%areaT(i,j), g%areaT(i+1,j))
463  if (abs(uh_col - sum(uh2d(i,:)))>uh_neglect) then
464  call mom_error(warning,"Column integral of uh does not match after "//&
465  "upwards redistribution")
466  endif
467 
468  enddo ! i-loop
469 
470  do k=1,nz ; do i=is-1,ie
471  uh(i,j,k) = uh2d(i,k)
472  enddo ; enddo
473  enddo
474 
475 end subroutine distribute_residual_uh_upwards
476 
477 !> In the case where offline advection has failed to converge, redistribute the u-flux
478 !! into layers above
479 subroutine distribute_residual_vh_upwards(G, GV, hvol, vh)
480  type(ocean_grid_type), pointer :: g !< ocean grid structure
481  type(verticalgrid_type), pointer :: gv !< ocean vertical grid structure
482  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
483  intent(in ) :: hvol !< Mass of water in the cells at the end
484  !! of the previous timestep [kg]
485  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
486  intent(inout) :: vh !< Meridional mass transport within a timestep [kg]
487 
488  real, dimension(SZJB_(G),SZK_(G)) :: vh2d
489  real, dimension(SZJB_(G)) :: vh2d_sum
490  real, dimension(SZJ_(G),SZK_(G)) :: h2d
491  real, dimension(SZJ_(G)) :: h2d_sum
492 
493  real :: vh_neglect, vh_remain, vh_col, vh_sum
494  real :: hup, hlos, min_h
495  integer :: i, j, k, m, is, ie, js, je, nz, k_rev
496 
497  ! Set index-related variables for fields on T-grid
498  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
499 
500  min_h = 0.1*gv%Angstrom_H
501 
502  do i=is,ie
503  ! Copy over uh and cell volume to working arrays
504  do k=1,nz ; do j=js-2,je+1
505  vh2d(j,k) = vh(i,j,k)
506  enddo ; enddo
507  do k=1,nz ; do j=js-1,je+1
508  h2d(j,k) = hvol(i,j,k)-min_h*g%US%L_to_m**2*g%areaT(i,j)
509  enddo ; enddo
510 
511  do j=js-1,je
512  vh_col = sum(vh2d(j,:))
513  do k=1,nz
514  vh_remain = vh2d(j,k)
515  vh2d(j,k) = 0.0
516  if (abs(vh_remain)>0.0) then
517  do k_rev = k,1,-1
518  vh_sum = vh_remain + vh2d(j,k_rev)
519  if (vh_sum<0.0) then ! Transport to the left
520  hup = h2d(j+1,k_rev)
521  hlos = max(0.0,vh2d(j+1,k_rev))
522  if ((((hup - hlos) + vh_sum) < 0.0) .and. &
523  ((0.5*hup + vh_sum) < 0.0)) then
524  vh2d(j,k_rev) = min(-0.5*hup,-hup+hlos,0.0)
525  vh_remain = vh_sum - vh2d(j,k_rev)
526  else
527  vh2d(j,k_rev) = vh_sum
528  vh_remain = 0.0
529  exit
530  endif
531  else ! Transport to the right
532  hup = h2d(j,k_rev)
533  hlos = max(0.0,-vh2d(j-1,k_rev))
534  if ((((hup - hlos) - vh_sum) < 0.0) .and. &
535  ((0.5*hup - vh_sum) < 0.0)) then
536  vh2d(j,k_rev) = max(0.5*hup,hup-hlos,0.0)
537  vh_remain = vh_sum - vh2d(j,k_rev)
538  else
539  vh2d(j,k_rev) = vh_sum
540  vh_remain = 0.0
541  exit
542  endif
543  endif
544 
545  enddo ! k_rev
546  endif
547 
548  if (abs(vh_remain)>0.0) then
549  if (k<nz) then
550  vh2d(j,k+1) = vh2d(j,k+1) + vh_remain
551  else
552  vh2d(j,k) = vh2d(j,k) + vh_remain
553  call mom_error(warning,"Water column cannot accommodate VH redistribution. Tracer will not be conserved")
554  endif
555  endif ! k-loop
556  enddo
557 
558  ! Calculate and check that column integrated transports match the original to
559  ! within the tolerance limit
560  vh_neglect = gv%Angstrom_H*g%US%L_to_m**2 * min(g%areaT(i,j), g%areaT(i,j+1))
561  if ( abs(vh_col-sum(vh2d(j,:))) > vh_neglect) then
562  call mom_error(warning,"Column integral of vh does not match after "//&
563  "upwards redistribution")
564  endif
565  enddo
566 
567  do k=1,nz ; do j=js-1,je
568  vh(i,j,k) = vh2d(j,k)
569  enddo ; enddo
570  enddo
571 
572 end subroutine distribute_residual_vh_upwards
573 
574 !> add_diurnal_SW adjusts the shortwave fluxes in an forcying_type variable
575 !! to add a synthetic diurnal cycle. Adapted from SIS2
576 subroutine offline_add_diurnal_sw(fluxes, G, Time_start, Time_end)
577  type(forcing), intent(inout) :: fluxes !< The type with atmospheric fluxes to be adjusted.
578  type(ocean_grid_type), intent(in) :: g !< The ocean lateral grid type.
579  type(time_type), intent(in) :: time_start !< The start time for this step.
580  type(time_type), intent(in) :: time_end !< The ending time for this step.
581 
582  real :: diurnal_factor, time_since_ae, rad
583  real :: fracday_dt, fracday_day
584  real :: cosz_day, cosz_dt, rrsun_day, rrsun_dt
585  type(time_type) :: dt_here
586 
587  integer :: i, j, k, i2, j2, isc, iec, jsc, jec, i_off, j_off
588 
589  isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
590  i_off = lbound(fluxes%sens,1) - g%isc ; j_off = lbound(fluxes%sens,2) - g%jsc
591 
592  ! Orbital_time extracts the time of year relative to the northern
593  ! hemisphere autumnal equinox from a time_type variable.
594  time_since_ae = orbital_time(time_start)
595  dt_here = time_end - time_start
596  rad = acos(-1.)/180.
597 
598 !$OMP parallel do default(none) shared(isc,iec,jsc,jec,G,rad,Time_start,dt_here,time_since_ae, &
599 !$OMP fluxes,i_off,j_off) &
600 !$OMP private(i,j,i2,j2,k,cosz_dt,fracday_dt,rrsun_dt, &
601 !$OMP fracday_day,cosz_day,rrsun_day,diurnal_factor)
602  do j=jsc,jec ; do i=isc,iec
603 ! Per Rick Hemler:
604 ! Call diurnal_solar with dtime=dt_here to get cosz averaged over dt_here.
605 ! Call daily_mean_solar to get cosz averaged over a day. Then
606 ! diurnal_factor = cosz_dt_ice*fracday_dt_ice*rrsun_dt_ice /
607 ! cosz_day*fracday_day*rrsun_day
608 
609  call diurnal_solar(g%geoLatT(i,j)*rad, g%geoLonT(i,j)*rad, time_start, cosz=cosz_dt, &
610  fracday=fracday_dt, rrsun=rrsun_dt, dt_time=dt_here)
611  call daily_mean_solar(g%geoLatT(i,j)*rad, time_since_ae, cosz_day, fracday_day, rrsun_day)
612  diurnal_factor = cosz_dt*fracday_dt*rrsun_dt / &
613  max(1e-30, cosz_day*fracday_day*rrsun_day)
614 
615  i2 = i+i_off ; j2 = j+j_off
616  fluxes%sw(i2,j2) = fluxes%sw(i2,j2) * diurnal_factor
617  fluxes%sw_vis_dir(i2,j2) = fluxes%sw_vis_dir(i2,j2) * diurnal_factor
618  fluxes%sw_vis_dif(i2,j2) = fluxes%sw_vis_dif(i2,j2) * diurnal_factor
619  fluxes%sw_nir_dir(i2,j2) = fluxes%sw_nir_dir(i2,j2) * diurnal_factor
620  fluxes%sw_nir_dif(i2,j2) = fluxes%sw_nir_dif(i2,j2) * diurnal_factor
621  enddo ; enddo
622 
623 end subroutine offline_add_diurnal_sw
624 
625 !> Controls the reading in 3d mass fluxes, diffusive fluxes, and other fields stored
626 !! in a previous integration of the online model
627 subroutine update_offline_from_files(G, GV, nk_input, mean_file, sum_file, snap_file, surf_file, h_end, &
628  uhtr, vhtr, temp_mean, salt_mean, mld, Kd, fluxes, ridx_sum, ridx_snap, read_mld, read_sw, &
629  read_ts_uvh, do_ale_in)
630 
631  type(ocean_grid_type), intent(inout) :: g !< Horizontal grid type
632  type(verticalgrid_type), intent(in ) :: gv !< Vertical grid type
633  integer, intent(in ) :: nk_input !< Number of levels in input file
634  character(len=*), intent(in ) :: mean_file !< Name of file with averages fields
635  character(len=*), intent(in ) :: sum_file !< Name of file with summed fields
636  character(len=*), intent(in ) :: snap_file !< Name of file with snapshot fields
637  character(len=*), intent(in ) :: surf_file !< Name of file with surface fields
638  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
639  intent(inout) :: uhtr !< Zonal mass fluxes [kg]
640  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
641  intent(inout) :: vhtr !< Meridional mass fluxes [kg]
642  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
643  intent(inout) :: h_end !< End of timestep layer thickness
644  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
645  intent(inout) :: temp_mean !< Averaged temperature
646  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
647  intent(inout) :: salt_mean !< Averaged salinity
648  real, dimension(SZI_(G),SZJ_(G)), &
649  intent(inout) :: mld !< Averaged mixed layer depth
650  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
651  intent(inout) :: kd !< Diapycnal diffusivities at interfaces
652  type(forcing), intent(inout) :: fluxes !< Fields with surface fluxes
653  integer, intent(in ) :: ridx_sum !< Read index for sum, mean, and surf files
654  integer, intent(in ) :: ridx_snap !< Read index for snapshot file
655  logical, intent(in ) :: read_mld !< True if reading in MLD
656  logical, intent(in ) :: read_sw !< True if reading in radiative fluxes
657  logical, intent(in ) :: read_ts_uvh !< True if reading in uh, vh, and h
658  logical, optional, intent(in ) :: do_ale_in !< True if using ALE algorithms
659 
660  logical :: do_ale
661  integer :: i, j, k, is, ie, js, je, nz
662  real :: initer_vert
663 
664  do_ale = .false.
665  if (present(do_ale_in) ) do_ale = do_ale_in
666 
667  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
668 
669  ! Check if reading in UH, VH, and h_end
670  if (read_ts_uvh) then
671  h_end(:,:,:) = 0.0
672  temp_mean(:,:,:) = 0.0
673  salt_mean(:,:,:) = 0.0
674  uhtr(:,:,:) = 0.0
675  vhtr(:,:,:) = 0.0
676  ! Time-summed fields
677  call mom_read_vector(sum_file, 'uhtr_sum', 'vhtr_sum', uhtr(:,:,1:nk_input), &
678  vhtr(:,:,1:nk_input), g%Domain, timelevel=ridx_sum)
679  call mom_read_data(snap_file, 'h_end', h_end(:,:,1:nk_input), g%Domain, &
680  timelevel=ridx_snap,position=center)
681  call mom_read_data(mean_file, 'temp', temp_mean(:,:,1:nk_input), g%Domain, &
682  timelevel=ridx_sum,position=center)
683  call mom_read_data(mean_file, 'salt', salt_mean(:,:,1:nk_input), g%Domain, &
684  timelevel=ridx_sum,position=center)
685  endif
686 
687  do j=js,je ; do i=is,ie
688  if (g%mask2dT(i,j)>0.) then
689  temp_mean(:,:,nk_input:nz) = temp_mean(i,j,nk_input)
690  salt_mean(:,:,nk_input:nz) = salt_mean(i,j,nk_input)
691  endif
692  enddo ; enddo
693 
694  ! Check if reading vertical diffusivities or entrainment fluxes
695  call mom_read_data( mean_file, 'Kd_interface', kd(:,:,1:nk_input+1), g%Domain, &
696  timelevel=ridx_sum,position=center)
697 
698  ! This block makes sure that the fluxes control structure, which may not be used in the solo_driver,
699  ! contains netMassIn and netMassOut which is necessary for the applyTracerBoundaryFluxesInOut routine
700  if (do_ale) then
701  if (.not. associated(fluxes%netMassOut)) then
702  allocate(fluxes%netMassOut(g%isd:g%ied,g%jsd:g%jed))
703  fluxes%netMassOut(:,:) = 0.0
704  endif
705  if (.not. associated(fluxes%netMassIn)) then
706  allocate(fluxes%netMassIn(g%isd:g%ied,g%jsd:g%jed))
707  fluxes%netMassIn(:,:) = 0.0
708  endif
709 
710  fluxes%netMassOut(:,:) = 0.0
711  fluxes%netMassIn(:,:) = 0.0
712  call mom_read_data(surf_file,'massout_flux_sum',fluxes%netMassOut, g%Domain, &
713  timelevel=ridx_sum)
714  call mom_read_data(surf_file,'massin_flux_sum', fluxes%netMassIn, g%Domain, &
715  timelevel=ridx_sum)
716 
717  do j=js,je ; do i=is,ie
718  if (g%mask2dT(i,j)<1.0) then
719  fluxes%netMassOut(i,j) = 0.0
720  fluxes%netMassIn(i,j) = 0.0
721  endif
722  enddo ; enddo
723 
724  endif
725 
726  if (read_mld) then
727  call mom_read_data(surf_file, 'ePBL_h_ML', mld, g%Domain, timelevel=ridx_sum)
728  endif
729 
730  if (read_sw) then
731  ! Shortwave radiation is only needed for offline mode with biogeochemistry but without the coupler.
732  ! Need to double check, but set_opacity seems to only need the sum of the diffuse and
733  ! direct fluxes in the visible and near-infrared bands. For convenience, we store the
734  ! sum of the direct and diffuse fluxes in the 'dir' field and set the 'dif' fields to zero
735  call mom_read_data(mean_file,'sw_vis', fluxes%sw_vis_dir, g%Domain, &
736  timelevel=ridx_sum, scale=g%US%W_m2_to_QRZ_T)
737  call mom_read_data(mean_file,'sw_nir', fluxes%sw_nir_dir, g%Domain, &
738  timelevel=ridx_sum, scale=g%US%W_m2_to_QRZ_T)
739  fluxes%sw_vis_dir(:,:) = fluxes%sw_vis_dir(:,:)*0.5
740  fluxes%sw_vis_dif(:,:) = fluxes%sw_vis_dir(:,:)
741  fluxes%sw_nir_dir(:,:) = fluxes%sw_nir_dir(:,:)*0.5
742  fluxes%sw_nir_dif(:,:) = fluxes%sw_nir_dir(:,:)
743  fluxes%sw = (fluxes%sw_vis_dir + fluxes%sw_vis_dif) + (fluxes%sw_nir_dir + fluxes%sw_nir_dif)
744  do j=js,je ; do i=is,ie
745  if (g%mask2dT(i,j)<1.0) then
746  fluxes%sw(i,j) = 0.0
747  fluxes%sw_vis_dir(i,j) = 0.0
748  fluxes%sw_nir_dir(i,j) = 0.0
749  fluxes%sw_vis_dif(i,j) = 0.0
750  fluxes%sw_nir_dif(i,j) = 0.0
751  endif
752  enddo ; enddo
753  call pass_var(fluxes%sw,g%Domain)
754  call pass_var(fluxes%sw_vis_dir,g%Domain)
755  call pass_var(fluxes%sw_vis_dif,g%Domain)
756  call pass_var(fluxes%sw_nir_dir,g%Domain)
757  call pass_var(fluxes%sw_nir_dif,g%Domain)
758  endif
759 
760 end subroutine update_offline_from_files
761 
762 !> Fields for offline transport are copied from the stored arrays read during initialization
763 subroutine update_offline_from_arrays(G, GV, nk_input, ridx_sum, mean_file, sum_file, snap_file, uhtr, vhtr, &
764  hend, uhtr_all, vhtr_all, hend_all, temp, salt, temp_all, salt_all )
765  type(ocean_grid_type), intent(inout) :: g !< Horizontal grid type
766  type(verticalgrid_type), intent(in ) :: gv !< Vertical grid type
767  integer, intent(in ) :: nk_input !< Number of levels in input file
768  integer, intent(in ) :: ridx_sum !< Index to read from
769  character(len=200), intent(in ) :: mean_file !< Name of file with averages fields
770  character(len=200), intent(in ) :: sum_file !< Name of file with summed fields
771  character(len=200), intent(in ) :: snap_file !< Name of file with snapshot fields
772  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(inout) :: uhtr !< Zonal mass fluxes [kg]
773  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(inout) :: vhtr !< Meridional mass fluxes [kg]
774  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: hend !< End of timestep layer thickness [kg m-2]
775  real, dimension(:,:,:,:), allocatable, intent(inout) :: uhtr_all !< Zonal mass fluxes [kg]
776  real, dimension(:,:,:,:), allocatable, intent(inout) :: vhtr_all !< Meridional mass fluxes [kg]
777  real, dimension(:,:,:,:), allocatable, intent(inout) :: hend_all !< End of timestep layer thickness [kg m-2]
778  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: temp !< Temperature array
779  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: salt !< Salinity array
780  real, dimension(:,:,:,:), allocatable, intent(inout) :: temp_all !< Temperature array
781  real, dimension(:,:,:,:), allocatable, intent(inout) :: salt_all !< Salinity array
782 
783  integer :: i, j, k, is, ie, js, je, nz
784  real, parameter :: fill_value = 0.
785  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
786 
787  ! Check that all fields are allocated (this is a redundant check)
788  if (.not. allocated(uhtr_all)) &
789  call mom_error(fatal, "uhtr_all not allocated before call to update_transport_from_arrays")
790  if (.not. allocated(vhtr_all)) &
791  call mom_error(fatal, "vhtr_all not allocated before call to update_transport_from_arrays")
792  if (.not. allocated(hend_all)) &
793  call mom_error(fatal, "hend_all not allocated before call to update_transport_from_arrays")
794  if (.not. allocated(temp_all)) &
795  call mom_error(fatal, "temp_all not allocated before call to update_transport_from_arrays")
796  if (.not. allocated(salt_all)) &
797  call mom_error(fatal, "salt_all not allocated before call to update_transport_from_arrays")
798 
799  ! Copy uh, vh, h_end, temp, and salt
800  do k=1,nk_input ; do j=js,je ; do i=is,ie
801  uhtr(i,j,k) = uhtr_all(i,j,k,ridx_sum)
802  vhtr(i,j,k) = vhtr_all(i,j,k,ridx_sum)
803  hend(i,j,k) = hend_all(i,j,k,ridx_sum)
804  temp(i,j,k) = temp_all(i,j,k,ridx_sum)
805  salt(i,j,k) = salt_all(i,j,k,ridx_sum)
806  enddo ; enddo ; enddo
807 
808  ! Fill the rest of the arrays with 0s (fill_value could probably be changed to a runtime parameter)
809  do k=nk_input+1,nz ; do j=js,je ; do i=is,ie
810  uhtr(i,j,k) = fill_value
811  vhtr(i,j,k) = fill_value
812  hend(i,j,k) = fill_value
813  temp(i,j,k) = fill_value
814  salt(i,j,k) = fill_value
815  enddo ; enddo ; enddo
816 
817 end subroutine update_offline_from_arrays
818 
819 !> Calculates the next timelevel to read from the input fields. This allows the 'looping'
820 !! of the fields
821 function next_modulo_time(inidx, numtime)
822  ! Returns the next time interval to be read
823  integer :: numtime ! Number of time levels in input fields
824  integer :: inidx ! The current time index
825 
826  integer :: read_index ! The index in the input files that corresponds
827  ! to the current timestep
828 
829  integer :: next_modulo_time
830 
831  read_index = mod(inidx+1,numtime)
832  if (read_index < 0) read_index = inidx-read_index
833  if (read_index == 0) read_index = numtime
834 
835  next_modulo_time = read_index
836 
837 end function next_modulo_time
838 
839 end module mom_offline_aux
840 
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_diag_vkernels
Provides kernels for single-column interpolation, re-integration (re-mapping of integrated quantities...
Definition: MOM_diag_vkernels.F90:3
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:54
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:59
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_variables::vertvisc_type
Vertical viscosities, drag coefficients, and related fields.
Definition: MOM_variables.F90:208
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_io::mom_read_vector
Read a pair of data fields representing the two components of a vector from a file.
Definition: MOM_io.F90:82
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_opacity::optics_type
This type is used to store information about ocean optical properties.
Definition: MOM_opacity.F90:25
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_offline_aux
Contains routines related to offline transport of tracers. These routines are likely to be called fro...
Definition: MOM_offline_aux.F90:3
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:66
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_opacity
Routines used to calculate the opacity of the ocean.
Definition: MOM_opacity.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26