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