MOM6
advection_test_tracer.F90
Go to the documentation of this file.
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
22 !********+*********+*********+*********+*********+*********+*********+**
23 !* *
24 !* By Robert Hallberg, 2002 *
25 !* *
26 !* This file contains an example of the code that is needed to set *
27 !* up and use a set (in this case eleven) of dynamically passive *
28 !* tracers. These tracers dye the inflowing water or water initially *
29 !* within a range of latitudes or water initially in a range of *
30 !* depths. *
31 !* *
32 !* A single subroutine is called from within each file to register *
33 !* each of the tracers for reinitialization and advection and to *
34 !* register the subroutine that initializes the tracers and set up *
35 !* their output and the subroutine that does any tracer physics or *
36 !* chemistry along with diapycnal mixing (included here because some *
37 !* tracers may float or swim vertically or dye diapycnal processes). *
38 !* *
39 !* *
40 !* Macros written all in capital letters are defined in MOM_memory.h. *
41 !* *
42 !* A small fragment of the grid is shown below: *
43 !* *
44 !* j+1 x ^ x ^ x At x: q *
45 !* j+1 > o > o > At ^: v *
46 !* j x ^ x ^ x At >: u *
47 !* j > o > o > At o: h, tr *
48 !* j-1 x ^ x ^ x *
49 !* i-1 i i+1 At x & ^: *
50 !* i i+1 At > & o: *
51 !* *
52 !* The boundaries always run through q grid points (x). *
53 !* *
54 !********+*********+*********+*********+*********+*********+*********+**
55 
56 use mom_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
57 use mom_diag_mediator, only : diag_ctrl
59 use mom_error_handler, only : mom_error, fatal, warning
61 use mom_forcing_type, only : forcing
62 use mom_grid, only : ocean_grid_type
63 use mom_hor_index, only : hor_index_type
64 use mom_io, only : file_exists, read_data, slasher, vardesc, var_desc, query_vardesc
68 use mom_time_manager, only : time_type, get_time
72 use mom_variables, only : surface
74 
75 use coupler_util, only : set_coupler_values, ind_csurf
77 
78 implicit none ; private
79 
80 #include <MOM_memory.h>
81 
85 
86 ! ntr is the number of tracers in this module.
87 integer, parameter :: ntr = 11
88 
89 type p3d
90  real, dimension(:,:,:), pointer :: p => null()
91 end type p3d
92 
93 type, public :: advection_test_tracer_cs ; private
94  integer :: ntr = ntr ! Number of tracers in this module
95  logical :: coupled_tracers = .false. ! These tracers are not offered to the
96  ! coupler.
97  character(len=200) :: tracer_ic_file ! The full path to the IC file, or " "
98  ! to initialize internally.
99  type(time_type), pointer :: time ! A pointer to the ocean model's clock.
100  type(tracer_registry_type), pointer :: tr_reg => null()
101  real, pointer :: tr(:,:,:,:) => null() ! The array of tracers used in this
102  ! subroutine, in g m-3?
103  real, pointer :: tr_aux(:,:,:,:) => null() ! The masked tracer concentration
104  ! for output, in g m-3.
105  type(p3d), dimension(NTR) :: &
106  tr_adx, &! Tracer zonal advective fluxes in g m-3 m3 s-1.
107  tr_ady, &! Tracer meridional advective fluxes in g m-3 m3 s-1.
108  tr_dfx, &! Tracer zonal diffusive fluxes in g m-3 m3 s-1.
109  tr_dfy ! Tracer meridional diffusive fluxes in g m-3 m3 s-1.
110  real :: land_val(ntr) = -1.0 ! The value of tr used where land is masked out.
111  logical :: mask_tracers ! If true, tracers are masked out in massless layers.
112  logical :: use_sponge
113  logical :: tracers_may_reinit
114 
115  real :: x_origin, x_width ! Parameters describing the test functions
116  real :: y_origin, y_width ! Parameters describing the test functions
117 
118  integer, dimension(NTR) :: ind_tr ! Indices returned by aof_set_coupler_flux
119  ! if it is used and the surface tracer concentrations are to be
120  ! provided to the coupler.
121 
122  type(diag_ctrl), pointer :: diag ! A structure that is used to regulate the
123  ! timing of diagnostic output.
124  type(mom_restart_cs), pointer :: restart_csp => null()
125 
126  integer, dimension(NTR) :: id_tracer = -1, id_tr_adx = -1, id_tr_ady = -1
127  integer, dimension(NTR) :: id_tr_dfx = -1, id_tr_dfy = -1
128 
129  type(vardesc) :: tr_desc(ntr)
131 
132 contains
133 
134 function register_advection_test_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
135  type(hor_index_type), intent(in) :: HI
136  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
137  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
138  type(advection_test_tracer_cs), pointer :: CS
139  type(tracer_registry_type), pointer :: tr_Reg
140  type(mom_restart_cs), pointer :: restart_CS
141 ! This subroutine is used to register tracer fields and subroutines
142 ! to be used with MOM.
143 ! Arguments: HI - A horizontal index type structure.
144 ! (in) GV - The ocean's vertical grid structure.
145 ! (in) param_file - A structure indicating the open file to parse for
146 ! model parameter values.
147 ! (in/out) CS - A pointer that is set to point to the control structure
148 ! for this module
149 ! (in/out) tr_Reg - A pointer that is set to point to the control structure
150 ! for the tracer advection and diffusion module.
151 ! (in) restart_CS - A pointer to the restart control structure.
152  character(len=80) :: name, longname
153 ! This include declares and sets the variable "version".
154 #include "version_variable.h"
155  character(len=40) :: mdl = "advection_test_tracer" ! This module's name.
156  character(len=200) :: inputdir
157  real, pointer :: tr_ptr(:,:,:) => null()
158  logical :: register_advection_test_tracer
159  integer :: isd, ied, jsd, jed, nz, m
160  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
161 
162  if (associated(cs)) then
163  call mom_error(warning, "register_advection_test_tracer called with an "// &
164  "associated control structure.")
165  return
166  endif
167  allocate(cs)
168 
169  ! Read all relevant parameters and write them to the model log.
170  call log_version(param_file, mdl, version, "")
171 
172  call get_param(param_file, mdl, "ADVECTION_TEST_X_ORIGIN", cs%x_origin, &
173  "The x-coorindate of the center of the test-functions.\n", default=0.)
174  call get_param(param_file, mdl, "ADVECTION_TEST_Y_ORIGIN", cs%y_origin, &
175  "The y-coorindate of the center of the test-functions.\n", default=0.)
176  call get_param(param_file, mdl, "ADVECTION_TEST_X_WIDTH", cs%x_width, &
177  "The x-width of the test-functions.\n", default=0.)
178  call get_param(param_file, mdl, "ADVECTION_TEST_Y_WIDTH", cs%y_width, &
179  "The y-width of the test-functions.\n", default=0.)
180  call get_param(param_file, mdl, "ADVECTION_TEST_TRACER_IC_FILE", cs%tracer_IC_file, &
181  "The name of a file from which to read the initial \n"//&
182  "conditions for the tracers, or blank to initialize \n"//&
183  "them internally.", default=" ")
184 
185  if (len_trim(cs%tracer_IC_file) >= 1) then
186  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
187  cs%tracer_IC_file = trim(slasher(inputdir))//trim(cs%tracer_IC_file)
188  call log_param(param_file, mdl, "INPUTDIR/ADVECTION_TEST_TRACER_IC_FILE", &
189  cs%tracer_IC_file)
190  endif
191  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
192  "If true, sponges may be applied anywhere in the domain. \n"//&
193  "The exact location and properties of those sponges are \n"//&
194  "specified from MOM_initialization.F90.", default=.false.)
195 
196  call get_param(param_file, mdl, "MASK_TRACERS_IN_MASSLESS_LAYERS", cs%mask_tracers, &
197  "If true, tracers will be masked out in massless layers. \n", &
198  default=.false.)
199  call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
200  "If true, tracers may go through the initialization code \n"//&
201  "if they are not found in the restart files. Otherwise \n"//&
202  "it is a fatal error if the tracers are not found in the \n"//&
203  "restart files of a restarted run.", default=.false.)
204 
205 
206  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
207  if (cs%mask_tracers) then
208  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
209  endif
210 
211  do m=1,ntr
212  if (m < 10) then ; write(name,'("tr",I1.1)') m
213  else ; write(name,'("tr",I2.2)') m ; endif
214  write(longname,'("Concentration of Tracer ",I2.2)') m
215  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
216 
217  ! This is needed to force the compiler not to do a copy in the registration
218  ! calls. Curses on the designers and implementers of Fortran90.
219  tr_ptr => cs%tr(:,:,:,m)
220  ! Register the tracer for the restart file.
221  call register_restart_field(tr_ptr, cs%tr_desc(m), &
222  .not. cs%tracers_may_reinit, restart_cs)
223  ! Register the tracer for horizontal advection & diffusion.
224  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
225  tr_desc_ptr=cs%tr_desc(m))
226 
227  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
228  ! values to the coupler (if any). This is meta-code and its arguments will
229  ! currently (deliberately) give fatal errors if it is used.
230  if (cs%coupled_tracers) &
231  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
232  flux_type=' ', implementation=' ', caller="register_advection_test_tracer")
233  enddo
234 
235  cs%tr_Reg => tr_reg
236  cs%restart_CSp => restart_cs
237  register_advection_test_tracer = .true.
239 
240 subroutine initialize_advection_test_tracer(restart, day, G, GV, h,diag, OBC, CS, &
241  sponge_CSp, diag_to_Z_CSp)
242  logical, intent(in) :: restart
243  type(time_type), target, intent(in) :: day
244  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
245  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
246  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
247  type(diag_ctrl), target, intent(in) :: diag
248  type(ocean_obc_type), pointer :: OBC
249  type(advection_test_tracer_cs), pointer :: CS
250  type(sponge_cs), pointer :: sponge_CSp
251  type(diag_to_z_cs), pointer :: diag_to_Z_CSp
252 ! This subroutine initializes the NTR tracer fields in tr(:,:,:,:)
253 ! and it sets up the tracer output.
254 
255 ! Arguments: restart - .true. if the fields have already been read from
256 ! a restart file.
257 ! (in) day - Time of the start of the run.
258 ! (in) G - The ocean's grid structure.
259 ! (in) GV - The ocean's vertical grid structure.
260 ! (in) h - Layer thickness, in m or kg m-2.
261 ! (in) diag - A structure that is used to regulate diagnostic output.
262 ! (in) OBC - This open boundary condition type specifies whether, where,
263 ! and what open boundary conditions are used.
264 ! (in/out) CS - The control structure returned by a previous call to
265 ! register_advection_test_tracer.
266 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
267 ! they are in use. Otherwise this may be unassociated.
268 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
269 ! in depth space.
270  real, allocatable :: temp(:,:,:)
271  real, pointer, dimension(:,:,:) :: &
272  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
273  obc_tr1_v => null() ! specify the values of tracer 1 that should come
274  ! in through u- and v- points through the open
275  ! boundary conditions, in the same units as tr.
276  character(len=16) :: name ! A variable's name in a NetCDF file.
277  character(len=72) :: longname ! The long name of that variable.
278  character(len=48) :: units ! The dimensions of the variable.
279  character(len=48) :: flux_units ! The units for tracer fluxes, usually
280  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
281  real, pointer :: tr_ptr(:,:,:) => null()
282  real :: PI ! 3.1415926... calculated as 4*atan(1)
283  real :: tr_y ! Initial zonally uniform tracer concentrations.
284  real :: dist2 ! The distance squared from a line, in m2.
285  real :: h_neglect ! A thickness that is so small it is usually lost
286  ! in roundoff and can be neglected, in m.
287  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
288  integer :: IsdB, IedB, JsdB, JedB
289  real :: tmpx, tmpy, locx, locy
290 
291  if (.not.associated(cs)) return
292  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
293  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
294  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
295  h_neglect = gv%H_subroundoff
296 
297  cs%diag => diag
298  cs%ntr = ntr
299  do m=1,ntr
300  call query_vardesc(cs%tr_desc(m), name=name, &
301  caller="initialize_advection_test_tracer")
302  if ((.not.restart) .or. (cs%tracers_may_reinit .and. .not. &
303  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
304  do k=1,nz ; do j=js,je ; do i=is,ie
305  cs%tr(i,j,k,m) = 0.0
306  enddo ; enddo ; enddo
307  k=1 ! Square wave
308  do j=js,je ; do i=is,ie
309  if (abs(g%geoLonT(i,j)-cs%x_origin)<0.5*cs%x_width .and. &
310  abs(g%geoLatT(i,j)-cs%y_origin)<0.5*cs%y_width) cs%tr(i,j,k,m) = 1.0
311  enddo ; enddo
312  k=2 ! Triangle wave
313  do j=js,je ; do i=is,ie
314  locx=abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
315  locy=abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
316  cs%tr(i,j,k,m) = max(0.0, 1.0-locx)*max(0.0, 1.0-locy)
317  enddo ; enddo
318  k=3 ! Cosine bell
319  do j=js,je ; do i=is,ie
320  locx=min(1.0, abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width)*(acos(0.0)*2.)
321  locy=min(1.0, abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width)*(acos(0.0)*2.)
322  cs%tr(i,j,k,m) = (1.0+cos(locx))*(1.0+cos(locy))*0.25
323  enddo ; enddo
324  k=4 ! Cylinder
325  do j=js,je ; do i=is,ie
326  locx=abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
327  locy=abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
328  if (locx**2+locy**2<=1.0) cs%tr(i,j,k,m) = 1.0
329  enddo ; enddo
330  k=5 ! Cut cylinder
331  do j=js,je ; do i=is,ie
332  locx=(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
333  locy=(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
334  if (locx**2+locy**2<=1.0) cs%tr(i,j,k,m) = 1.0
335  if (locx>0.0.and.abs(locy)<0.2) cs%tr(i,j,k,m) = 0.0
336  enddo ; enddo
337  endif ! restart
338  enddo
339 
340 
341  ! This needs to be changed if the units of tracer are changed above.
342  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
343  else ; flux_units = "kg s-1" ; endif
344 
345  do m=1,ntr
346  ! Register the tracer for the restart file.
347  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
348  caller="initialize_advection_test_tracer")
349  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
350  day, trim(longname) , trim(units))
351  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
352  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
353  trim(flux_units))
354  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
355  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
356  trim(flux_units))
357  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
358  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
359  trim(flux_units))
360  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
361  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
362  trim(flux_units))
363  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
364  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
365  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
366  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
367 
368 ! Register the tracer for horizontal advection & diffusion.
369  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
370  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
371  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
372  cs%tr_ady(m)%p,cs%tr_dfx(m)%p,cs%tr_dfy(m)%p)
373 
374  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
375  day, g, diag_to_z_csp)
376  enddo
377 
379 
380 
381 subroutine advection_test_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
382  evap_CFL_limit, minimum_forcing_depth)
383  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
384  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
385  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
386  type(forcing), intent(in) :: fluxes
387  real, intent(in) :: dt !< The amount of time covered by this call, in s
388  type(advection_test_tracer_cs), pointer :: CS
389  real, optional,intent(in) :: evap_CFL_limit
390  real, optional,intent(in) :: minimum_forcing_depth
391 ! This subroutine applies diapycnal diffusion and any other column
392 ! tracer physics or chemistry to the tracers from this file.
393 ! This is a simple example of a set of advected passive tracers.
394 
395 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
396 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
397 ! (in) ea - an array to which the amount of fluid entrained
398 ! from the layer above during this call will be
399 ! added, in m or kg m-2.
400 ! (in) eb - an array to which the amount of fluid entrained
401 ! from the layer below during this call will be
402 ! added, in m or kg m-2.
403 ! (in) fluxes - A structure containing pointers to any possible
404 ! forcing fields. Unused fields have NULL ptrs.
405 ! (in) dt - The amount of time covered by this call, in s.
406 ! (in) G - The ocean's grid structure.
407 ! (in) GV - The ocean's vertical grid structure.
408 ! (in) CS - The control structure returned by a previous call to
409 ! register_advection_test_tracer.
410 !
411 ! The arguments to this subroutine are redundant in that
412 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
413  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
414  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
415  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
416  integer :: i, j, k, is, ie, js, je, nz, m
417  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
418 
419  if (.not.associated(cs)) return
420 
421  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
422  do m=1,cs%ntr
423  do k=1,nz ;do j=js,je ; do i=is,ie
424  h_work(i,j,k) = h_old(i,j,k)
425  enddo ; enddo ; enddo;
426  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
427  evap_cfl_limit, minimum_forcing_depth)
428  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
429  enddo
430  else
431  do m=1,ntr
432  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
433  enddo
434  endif
435 
436  if (cs%mask_tracers) then
437  do m = 1,ntr ; if (cs%id_tracer(m) > 0) then
438  do k=1,nz ; do j=js,je ; do i=is,ie
439  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
440  cs%tr_aux(i,j,k,m) = cs%land_val(m)
441  else
442  cs%tr_aux(i,j,k,m) = cs%tr(i,j,k,m)
443  endif
444  enddo ; enddo ; enddo
445  endif ; enddo
446  endif
447 
448  do m=1,ntr
449  if (cs%mask_tracers) then
450  if (cs%id_tracer(m)>0) &
451  call post_data(cs%id_tracer(m),cs%tr_aux(:,:,:,m),cs%diag)
452  else
453  if (cs%id_tracer(m)>0) &
454  call post_data(cs%id_tracer(m),cs%tr(:,:,:,m),cs%diag)
455  endif
456  if (cs%id_tr_adx(m)>0) &
457  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
458  if (cs%id_tr_ady(m)>0) &
459  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
460  if (cs%id_tr_dfx(m)>0) &
461  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
462  if (cs%id_tr_dfy(m)>0) &
463  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
464  enddo
466 
467 subroutine advection_test_tracer_surface_state(state, h, G, CS)
468  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
469  type(surface), intent(inout) :: state
470  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
471  type(advection_test_tracer_cs), pointer :: CS
472 ! This particular tracer package does not report anything back to the coupler.
473 ! The code that is here is just a rough guide for packages that would.
474 ! Arguments: state - A structure containing fields that describe the
475 ! surface state of the ocean.
476 ! (in) h - Layer thickness, in m or kg m-2.
477 ! (in) G - The ocean's grid structure.
478 ! (in) CS - The control structure returned by a previous call to
479 ! register_advection_test_tracer.
480  integer :: i, j, m, is, ie, js, je
481  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
482 
483  if (.not.associated(cs)) return
484 
485  if (cs%coupled_tracers) then
486  do m=1,ntr
487  ! This call loads the surface vlues into the appropriate array in the
488  ! coupler-type structure.
489  call set_coupler_values(cs%tr(:,:,1,1), state%tr_fields, cs%ind_tr(m), &
490  ind_csurf, is, ie, js, je)
491  enddo
492  endif
494 
495 function advection_test_stock(h, stocks, G, GV, CS, names, units, stock_index)
496  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
497  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
498  real, dimension(:), intent(out) :: stocks
499  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
500  type(advection_test_tracer_cs), pointer :: CS
501  character(len=*), dimension(:), intent(out) :: names
502  character(len=*), dimension(:), intent(out) :: units
503  integer, optional, intent(in) :: stock_index
504  integer :: advection_test_stock
505 ! This function calculates the mass-weighted integral of all tracer stocks,
506 ! returning the number of stocks it has calculated. If the stock_index
507 ! is present, only the stock corresponding to that coded index is returned.
508 
509 ! Arguments: h - Layer thickness, in m or kg m-2.
510 ! (out) stocks - the mass-weighted integrated amount of each tracer,
511 ! in kg times concentration units.
512 ! (in) G - The ocean's grid structure.
513 ! (in) GV - The ocean's vertical grid structure.
514 ! (in) CS - The control structure returned by a previous call to
515 ! register_ideal_age_tracer.
516 ! (out) names - the names of the stocks calculated.
517 ! (out) units - the units of the stocks calculated.
518 ! (in,opt) stock_index - the coded index of a specific stock being sought.
519 ! Return value: the number of stocks calculated here.
520 
521  integer :: i, j, k, is, ie, js, je, nz, m
522  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
523 
524  advection_test_stock = 0
525  if (.not.associated(cs)) return
526  if (cs%ntr < 1) return
527 
528  if (present(stock_index)) then ; if (stock_index > 0) then
529  ! Check whether this stock is available from this routine.
530 
531  ! No stocks from this routine are being checked yet. Return 0.
532  return
533  endif ; endif
534 
535  do m=1,cs%ntr
536  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="advection_test_stock")
537  stocks(m) = 0.0
538  do k=1,nz ; do j=js,je ; do i=is,ie
539  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
540  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
541  enddo ; enddo ; enddo
542  stocks(m) = gv%H_to_kg_m2 * stocks(m)
543  enddo
544  advection_test_stock = cs%ntr
545 
546 end function advection_test_stock
547 
548 subroutine advection_test_tracer_end(CS)
549  type(advection_test_tracer_cs), pointer :: CS
550  integer :: m
551 
552  if (associated(cs)) then
553  if (associated(cs%tr)) deallocate(cs%tr)
554  if (associated(cs%tr_aux)) deallocate(cs%tr_aux)
555  do m=1,ntr
556  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
557  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
558  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
559  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
560  enddo
561 
562  deallocate(cs)
563  endif
564 end subroutine advection_test_tracer_end
565 
566 end module advection_test_tracer
The following structure contains pointers to various fields which may be used describe the surface st...
subroutine, public advection_test_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
Definition: MOM_io.F90:585
This module implements boundary forcing for MOM6.
integer function, public aof_set_coupler_flux(name, flux_type, implementation, atm_tr_index, param, flag, ice_restart_file, ocean_restart_file, units, caller)
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Provides the ocean grid type.
Definition: MOM_grid.F90:2
integer function, public advection_test_stock(h, stocks, G, GV, CS, names, units, stock_index)
This module contains I/O framework code.
Definition: MOM_io.F90:2
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Container for horizontal index ranges for data, computational and global domains. ...
subroutine, public register_tracer(tr1, tr_desc, param_file, HI, GV, Reg, tr_desc_ptr, ad_x, ad_y, df_x, df_y, OBC_inflow, OBC_in_u, OBC_in_v, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine registers a tracer to be advected and laterally diffused.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
subroutine, public applytracerboundaryfluxesinout(G, GV, Tr, dt, fluxes, h, evap_CFL_limit, minimum_forcing_depth, in_flux_optional, out_flux_optional, update_h_opt)
This routine is modeled after applyBoundaryFluxesInOut in MOM_diabatic_aux.F90 NOTE: Please note that...
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
subroutine, public register_z_tracer(tr_ptr, name, long_name, units, Time, G, CS, standard_name, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name)
This subroutine registers a tracer to be output in depth space.
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public add_tracer_diagnostics(name, Reg, ad_x, ad_y, df_x, df_y, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine adds diagnostic arrays for a tracer that has previously been registered by a call to ...
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public set_coupler_values(array_in, BC_struc, BC_index, BC_element, is, ie, js, je, conversion)
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
subroutine, public tracer_vertdiff(h_old, ea, eb, dt, tr, G, GV, sfc_flux, btm_flux, btm_reservoir, sink_rate, convert_flux_in)
This subroutine solves a tridiagonal equation for the final tracer concentrations after the dual-entr...
logical function, public register_advection_test_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
subroutine, public advection_test_tracer_end(CS)
subroutine, public query_vardesc(vd, name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
This routine queries vardesc.
Definition: MOM_io.F90:664
subroutine, public initialize_advection_test_tracer(restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp)
Controls where open boundary conditions are applied.
subroutine, public advection_test_tracer_surface_state(state, h, G, CS)
integer function, public register_diag_field(module_name, field_name, axes, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, cell_methods, x_cell_method, y_cell_method, v_cell_method, conversion, v_extensive)
Returns the "diag_mediator" handle for a group (native, CMOR, z-coord, ...) of diagnostics derived fr...
subroutine, public mom_error(level, message, all_print)