MOM6
tracer_example.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 one) of dynamically passive tracers.*
28 !* *
29 !* A single subroutine is called from within each file to register *
30 !* each of the tracers for reinitialization and advection and to *
31 !* register the subroutine that initializes the tracers and set up *
32 !* their output and the subroutine that does any tracer physics or *
33 !* chemistry along with diapycnal mixing (included here because some *
34 !* tracers may float or swim vertically or dye diapycnal processes). *
35 !* *
36 !* *
37 !* Macros written all in capital letters are defined in MOM_memory.h. *
38 !* *
39 !* A small fragment of the grid is shown below: *
40 !* *
41 !* j+1 x ^ x ^ x At x: q *
42 !* j+1 > o > o > At ^: v *
43 !* j x ^ x ^ x At >: u *
44 !* j > o > o > At o: h, tr *
45 !* j-1 x ^ x ^ x *
46 !* i-1 i i+1 At x & ^: *
47 !* i i+1 At > & o: *
48 !* *
49 !* The boundaries always run through q grid points (x). *
50 !* *
51 !********+*********+*********+*********+*********+*********+*********+**
52 
53 use mom_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
54 use mom_diag_mediator, only : diag_ctrl
56 use mom_error_handler, only : mom_error, fatal, warning
58 use mom_forcing_type, only : forcing
59 use mom_grid, only : ocean_grid_type
60 use mom_hor_index, only : hor_index_type
61 use mom_io, only : file_exists, read_data, slasher, vardesc, var_desc, query_vardesc
65 use mom_time_manager, only : time_type, get_time
68 use mom_variables, only : surface
70 
71 use coupler_util, only : set_coupler_values, ind_csurf
73 
74 implicit none ; private
75 
76 #include <MOM_memory.h>
77 
80 
81 ! NTR is the number of tracers in this module.
82 integer, parameter :: ntr = 1
83 
84 type p3d
85  real, dimension(:,:,:), pointer :: p => null()
86 end type p3d
87 
88 type, public :: user_tracer_example_cs ; private
89  logical :: coupled_tracers = .false. ! These tracers are not offered to the
90  ! coupler.
91  character(len=200) :: tracer_ic_file ! The full path to the IC file, or " "
92  ! to initialize internally.
93  type(time_type), pointer :: time ! A pointer to the ocean model's clock.
94  type(tracer_registry_type), pointer :: tr_reg => null()
95  real, pointer :: tr(:,:,:,:) => null() ! The array of tracers used in this
96  ! subroutine, in g m-3?
97  real, pointer :: tr_aux(:,:,:,:) => null() ! The masked tracer concentration
98  ! for output, in g m-3.
99  type(p3d), dimension(NTR) :: &
100  tr_adx, &! Tracer zonal advective fluxes in g m-3 m3 s-1.
101  tr_ady, &! Tracer meridional advective fluxes in g m-3 m3 s-1.
102  tr_dfx, &! Tracer zonal diffusive fluxes in g m-3 m3 s-1.
103  tr_dfy ! Tracer meridional diffusive fluxes in g m-3 m3 s-1.
104  real :: land_val(ntr) = -1.0 ! The value of tr used where land is masked out.
105  logical :: mask_tracers ! If true, tracers are masked out in massless layers.
106  logical :: use_sponge
107 
108  integer, dimension(NTR) :: ind_tr ! Indices returned by aof_set_coupler_flux
109  ! if it is used and the surface tracer concentrations are to be
110  ! provided to the coupler.
111 
112  type(diag_ctrl), pointer :: diag ! A pointer to a structure of shareable
113  ! ocean diagnostic fields and control variables.
114  integer, dimension(NTR) :: id_tracer = -1, id_tr_adx = -1, id_tr_ady = -1
115  integer, dimension(NTR) :: id_tr_dfx = -1, id_tr_dfy = -1
116 
117  type(vardesc) :: tr_desc(ntr)
118 end type user_tracer_example_cs
119 
120 contains
121 
122 function user_register_tracer_example(HI, GV, param_file, CS, tr_Reg, restart_CS)
123  type(hor_index_type), intent(in) :: HI
124  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
125  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
126  type(user_tracer_example_cs), pointer :: CS
127  type(tracer_registry_type), pointer :: tr_Reg
128  type(mom_restart_cs), pointer :: restart_CS
129 ! This subroutine is used to register tracer fields and subroutines
130 ! to be used with MOM.
131 ! Arguments: HI - A horizontal index type structure.
132 ! (in) GV - The ocean's vertical grid structure.
133 ! (in) param_file - A structure indicating the open file to parse for
134 ! model parameter values.
135 ! (in/out) CS - A pointer that is set to point to the control structure
136 ! for this module
137 ! (in/out) tr_Reg - A pointer that is set to point to the control structure
138 ! for the tracer advection and diffusion module.
139 ! (in) restart_CS - A pointer to the restart control structure.
140  character(len=80) :: name, longname
141 ! This include declares and sets the variable "version".
142 #include "version_variable.h"
143  character(len=40) :: mdl = "tracer_example" ! This module's name.
144  character(len=200) :: inputdir
145  real, pointer :: tr_ptr(:,:,:) => null()
146  logical :: USER_register_tracer_example
147  integer :: isd, ied, jsd, jed, nz, m
148  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
149 
150  if (associated(cs)) then
151  call mom_error(warning, "USER_register_tracer_example called with an "// &
152  "associated control structure.")
153  return
154  endif
155  allocate(cs)
156 
157  ! Read all relevant parameters and write them to the model log.
158  call log_version(param_file, mdl, version, "")
159  call get_param(param_file, mdl, "TRACER_EXAMPLE_IC_FILE", cs%tracer_IC_file, &
160  "The name of a file from which to read the initial \n"//&
161  "conditions for the DOME tracers, or blank to initialize \n"//&
162  "them internally.", default=" ")
163  if (len_trim(cs%tracer_IC_file) >= 1) then
164  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
165  cs%tracer_IC_file = trim(slasher(inputdir))//trim(cs%tracer_IC_file)
166  call log_param(param_file, mdl, "INPUTDIR/TRACER_EXAMPLE_IC_FILE", &
167  cs%tracer_IC_file)
168  endif
169  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
170  "If true, sponges may be applied anywhere in the domain. \n"//&
171  "The exact location and properties of those sponges are \n"//&
172  "specified from MOM_initialization.F90.", default=.false.)
173 
174  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
175  if (cs%mask_tracers) then
176  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
177  endif
178 
179  do m=1,ntr
180  if (m < 10) then ; write(name,'("tr",I1.1)') m
181  else ; write(name,'("tr",I2.2)') m ; endif
182  write(longname,'("Concentration of Tracer ",I2.2)') m
183  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
184 
185  ! This is needed to force the compiler not to do a copy in the registration
186  ! calls. Curses on the designers and implementers of Fortran90.
187  tr_ptr => cs%tr(:,:,:,m)
188  ! Register the tracer for the restart file.
189  call register_restart_field(tr_ptr, cs%tr_desc(m), .true., restart_cs)
190  ! Register the tracer for horizontal advection & diffusion.
191  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
192  tr_desc_ptr=cs%tr_desc(m))
193 
194  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
195  ! values to the coupler (if any). This is meta-code and its arguments will
196  ! currently (deliberately) give fatal errors if it is used.
197  if (cs%coupled_tracers) &
198  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
199  flux_type=' ', implementation=' ', caller="USER_register_tracer_example")
200  enddo
201 
202  cs%tr_Reg => tr_reg
203  user_register_tracer_example = .true.
204 end function user_register_tracer_example
205 
206 subroutine user_initialize_tracer(restart, day, G, GV, h, diag, OBC, CS, &
207  sponge_CSp, diag_to_Z_CSp)
208  logical, intent(in) :: restart
209  type(time_type), target, intent(in) :: day
210  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
211  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
212  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
213  intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
214  type(diag_ctrl), target, intent(in) :: diag
215  type(ocean_obc_type), pointer :: OBC
216  type(user_tracer_example_cs), pointer :: CS
217  type(sponge_cs), pointer :: sponge_CSp
218  type(diag_to_z_cs), pointer :: diag_to_Z_CSp
219 ! This subroutine initializes the NTR tracer fields in tr(:,:,:,:)
220 ! and it sets up the tracer output.
221 
222 ! Arguments: restart - .true. if the fields have already been read from
223 ! a restart file.
224 ! (in) day - Time of the start of the run.
225 ! (in) G - The ocean's grid structure.
226 ! (in) GV - The ocean's vertical grid structure.
227 ! (in) h - Layer thickness, in m or kg m-2.
228 ! (in) diag - A structure that is used to regulate diagnostic output.
229 ! (in) OBC - This open boundary condition type specifies whether, where,
230 ! and what open boundary conditions are used.
231 ! (in/out) CS - The control structure returned by a previous call to
232 ! USER_register_tracer_example.
233 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
234 ! they are in use. Otherwise this may be unassociated.
235 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
236 ! in depth space.
237  real, allocatable :: temp(:,:,:)
238  real, pointer, dimension(:,:,:) :: &
239  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
240  obc_tr1_v => null() ! specify the values of tracer 1 that should come
241  ! in through u- and v- points through the open
242  ! boundary conditions, in the same units as tr.
243  character(len=32) :: name ! A variable's name in a NetCDF file.
244  character(len=72) :: longname ! The long name of that variable.
245  character(len=48) :: units ! The dimensions of the variable.
246  character(len=48) :: flux_units ! The units for tracer fluxes, usually
247  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
248  real, pointer :: tr_ptr(:,:,:) => null()
249  real :: PI ! 3.1415926... calculated as 4*atan(1)
250  real :: tr_y ! Initial zonally uniform tracer concentrations.
251  real :: dist2 ! The distance squared from a line, in m2.
252  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
253  integer :: IsdB, IedB, JsdB, JedB, lntr
254 
255  if (.not.associated(cs)) return
256  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
257  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
258  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
259 
260  lntr = ntr ! Avoids compile-time warning when NTR<2
261  cs%Time => day
262  cs%diag => diag
263 
264  if (.not.restart) then
265  if (len_trim(cs%tracer_IC_file) >= 1) then
266 ! Read the tracer concentrations from a netcdf file.
267  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
268  call mom_error(fatal, "USER_initialize_tracer: Unable to open "// &
269  cs%tracer_IC_file)
270  do m=1,ntr
271  call query_vardesc(cs%tr_desc(m), name, caller="USER_initialize_tracer")
272  call read_data(cs%tracer_IC_file, trim(name), &
273  cs%tr(:,:,:,m), domain=g%Domain%mpp_domain)
274  enddo
275  else
276  do m=1,ntr
277  do k=1,nz ; do j=js,je ; do i=is,ie
278  cs%tr(i,j,k,m) = 1.0e-20 ! This could just as well be 0.
279  enddo ; enddo ; enddo
280  enddo
281 
282 ! This sets a stripe of tracer across the basin.
283  pi = 4.0*atan(1.0)
284  do j=js,je
285  dist2 = (g%Rad_Earth * pi / 180.0)**2 * &
286  (g%geoLatT(i,j) - 40.0) * (g%geoLatT(i,j) - 40.0)
287  tr_y = 0.5*exp(-dist2/(1.0e5*1.0e5))
288 
289  do k=1,nz ; do i=is,ie
290 ! This adds the stripes of tracer to every layer.
291  cs%tr(i,j,k,1) = cs%tr(i,j,k,1) + tr_y
292  enddo ; enddo
293  enddo
294  endif
295  endif ! restart
296 
297  if ( cs%use_sponge ) then
298 ! If sponges are used, this example damps tracers in sponges in the
299 ! northern half of the domain to 1 and tracers in the southern half
300 ! to 0. For any tracers that are not damped in the sponge, the call
301 ! to set_up_sponge_field can simply be omitted.
302  if (.not.associated(sponge_csp)) &
303  call mom_error(fatal, "USER_initialize_tracer: "// &
304  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
305 
306  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
307  do k=1,nz ; do j=js,je ; do i=is,ie
308  if (g%geoLatT(i,j) > 700.0 .and. (k > nz/2)) then
309  temp(i,j,k) = 1.0
310  else
311  temp(i,j,k) = 0.0
312  endif
313  enddo ; enddo ; enddo
314 
315 ! do m=1,NTR
316  do m=1,1
317  ! This is needed to force the compiler not to do a copy in the sponge
318  ! calls. Curses on the designers and implementers of Fortran90.
319  tr_ptr => cs%tr(:,:,:,m)
320  call set_up_sponge_field(temp, tr_ptr, g, nz, sponge_csp)
321  enddo
322  deallocate(temp)
323  endif
324 
325  if (associated(obc)) then
326  call query_vardesc(cs%tr_desc(1), name, caller="USER_initialize_tracer")
327  if (obc%specified_v_BCs_exist_globally) then
328  allocate(obc_tr1_v(g%isd:g%ied,g%jsd:g%jed,nz))
329  do k=1,nz ; do j=g%jsd,g%jed ; do i=g%isd,g%ied
330  if (k < nz/2) then ; obc_tr1_v(i,j,k) = 0.0
331  else ; obc_tr1_v(i,j,k) = 1.0 ; endif
332  enddo ; enddo ; enddo
333  call add_tracer_obc_values(trim(name), cs%tr_Reg, &
334  0.0, obc_in_v=obc_tr1_v)
335  else
336  ! This is not expected in the DOME example.
337  call add_tracer_obc_values(trim(name), cs%tr_Reg, 0.0)
338  endif
339  ! All tracers but the first have 0 concentration in their inflows. As this
340  ! is the default value, the following calls are unnecessary.
341  do m=2,lntr
342  call query_vardesc(cs%tr_desc(m), name, caller="USER_initialize_tracer")
343  call add_tracer_obc_values(trim(name), cs%tr_Reg, 0.0)
344  enddo
345  endif
346 
347  ! This needs to be changed if the units of tracer are changed above.
348  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
349  else ; flux_units = "kg s-1" ; endif
350 
351  do m=1,ntr
352  ! Register the tracer for the restart file.
353  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
354  caller="USER_initialize_tracer")
355  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
356  day, trim(longname) , trim(units))
357  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
358  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
359  trim(flux_units))
360  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
361  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
362  trim(flux_units))
363  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
364  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
365  trim(flux_units))
366  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
367  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
368  trim(flux_units))
369  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
370  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
371  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
372  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
373 
374 ! Register the tracer for horizontal advection & diffusion.
375  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
376  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
377  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
378  cs%tr_ady(m)%p,cs%tr_dfx(m)%p,cs%tr_dfy(m)%p)
379 
380  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
381  day, g, diag_to_z_csp)
382  enddo
383 
384 end subroutine user_initialize_tracer
385 
386 subroutine tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS)
387  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
388  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
389  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
390  type(forcing), intent(in) :: fluxes
391  real, intent(in) :: dt !< The amount of time covered by this call, in s
392  type(user_tracer_example_cs), pointer :: CS
393 ! This subroutine applies diapycnal diffusion and any other column
394 ! tracer physics or chemistry to the tracers from this file.
395 ! This is a simple example of a set of advected passive tracers.
396 
397 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
398 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
399 ! (in) ea - an array to which the amount of fluid entrained
400 ! from the layer above during this call will be
401 ! added, in m or kg m-2.
402 ! (in) eb - an array to which the amount of fluid entrained
403 ! from the layer below during this call will be
404 ! added, in m or kg m-2.
405 ! (in) fluxes - A structure containing pointers to any possible
406 ! forcing fields. Unused fields have NULL ptrs.
407 ! (in) dt - The amount of time covered by this call, in s.
408 ! (in) G - The ocean's grid structure.
409 ! (in) GV - The ocean's vertical grid structure.
410 ! (in) CS - The control structure returned by a previous call to
411 ! USER_register_tracer_example.
412 !
413 ! The arguments to this subroutine are redundant in that
414 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
415 
416  real :: hold0(szi_(g)) ! The original topmost layer thickness,
417  ! with surface mass fluxes added back, m.
418  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
419  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
420  real :: d1(szi_(g)) ! d1=1-c1 is used by the tridiagonal solver.
421  real :: h_neglect ! A thickness that is so small it is usually lost
422  ! in roundoff and can be neglected, in m.
423  real :: b_denom_1 ! The first term in the denominator of b1, in m or kg m-2.
424  integer :: i, j, k, is, ie, js, je, nz, m
425 
426 ! The following array (trdc) determines the behavior of the tracer
427 ! diapycnal advection. The first element is 1 if tracers are
428 ! passively advected. The second and third are the concentrations
429 ! to which downwelling and upwelling water are set, respectively.
430 ! For most (normal) tracers, the appropriate vales are {1,0,0}.
431 
432  real :: trdc(3)
433 ! Uncomment the following line to dye both upwelling and downwelling.
434 ! data trdc / 0.0,1.0,1.0 /
435 ! Uncomment the following line to dye downwelling.
436 ! data trdc / 0.0,1.0,0.0 /
437 ! Uncomment the following line to dye upwelling.
438 ! data trdc / 0.0,0.0,1.0 /
439 ! Uncomment the following line for tracer concentrations to be set
440 ! to zero in any diapycnal motions.
441 ! data trdc / 0.0,0.0,0.0 /
442 ! Uncomment the following line for most "physical" tracers, which
443 ! are advected diapycnally in the usual manner.
444  data trdc / 1.0,0.0,0.0 /
445  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
446 
447  if (.not.associated(cs)) return
448  h_neglect = gv%H_subroundoff
449 
450  do j=js,je
451  do i=is,ie
452 ! The following line is appropriate for quantities like salinity
453 ! that are left behind by evaporation, and any surface fluxes would
454 ! be explicitly included in the flux structure.
455  hold0(i) = h_old(i,j,1)
456 ! The following line is appropriate for quantities like temperature
457 ! that can be assumed to have the same concentration in evaporation
458 ! as they had in the water. The explicit surface fluxes here would
459 ! reflect differences in concentration from the ambient water, not
460 ! the absolute fluxes.
461  ! hold0(i) = h_old(i,j,1) + ea(i,j,1)
462  b_denom_1 = h_old(i,j,1) + ea(i,j,1) + h_neglect
463  b1(i) = 1.0 / (b_denom_1 + eb(i,j,1))
464 ! d1(i) = b_denom_1 * b1(i)
465  d1(i) = trdc(1) * (b_denom_1 * b1(i)) + (1.0 - trdc(1))
466  do m=1,ntr
467  cs%tr(i,j,1,m) = b1(i)*(hold0(i)*cs%tr(i,j,1,m) + trdc(3)*eb(i,j,1))
468  ! Add any surface tracer fluxes to the preceding line.
469  enddo
470  enddo
471  do k=2,nz ; do i=is,ie
472  c1(i,k) = trdc(1) * eb(i,j,k-1) * b1(i)
473  b_denom_1 = h_old(i,j,k) + d1(i)*ea(i,j,k) + h_neglect
474  b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
475  d1(i) = trdc(1) * (b_denom_1 * b1(i)) + (1.0 - trdc(1))
476  do m=1,ntr
477  cs%tr(i,j,k,m) = b1(i) * (h_old(i,j,k)*cs%tr(i,j,k,m) + &
478  ea(i,j,k)*(trdc(1)*cs%tr(i,j,k-1,m)+trdc(2)) + &
479  eb(i,j,k)*trdc(3))
480  enddo
481  enddo ; enddo
482  do m=1,ntr ; do k=nz-1,1,-1 ; do i=is,ie
483  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + c1(i,k+1)*cs%tr(i,j,k+1,m)
484  enddo ; enddo ; enddo
485  enddo
486 
487  if (cs%mask_tracers) then
488  do m = 1,ntr ; if (cs%id_tracer(m) > 0) then
489  do k=1,nz ; do j=js,je ; do i=is,ie
490  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
491  cs%tr_aux(i,j,k,m) = cs%land_val(m)
492  else
493  cs%tr_aux(i,j,k,m) = cs%tr(i,j,k,m)
494  endif
495  enddo ; enddo ; enddo
496  endif ; enddo
497  endif
498 
499  do m=1,ntr
500  if (cs%mask_tracers) then
501  if (cs%id_tracer(m)>0) &
502  call post_data(cs%id_tracer(m),cs%tr_aux(:,:,:,m),cs%diag)
503  else
504  if (cs%id_tracer(m)>0) &
505  call post_data(cs%id_tracer(m),cs%tr(:,:,:,m),cs%diag)
506  endif
507  if (cs%id_tr_adx(m)>0) &
508  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
509  if (cs%id_tr_ady(m)>0) &
510  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
511  if (cs%id_tr_dfx(m)>0) &
512  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
513  if (cs%id_tr_dfy(m)>0) &
514  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
515  enddo
516 
517 end subroutine tracer_column_physics
518 
519 function user_tracer_stock(h, stocks, G, GV, CS, names, units, stock_index)
520  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
521  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
522  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
523  intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
524  real, dimension(:), intent(out) :: stocks
525  type(user_tracer_example_cs), pointer :: CS
526  character(len=*), dimension(:), intent(out) :: names
527  character(len=*), dimension(:), intent(out) :: units
528  integer, optional, intent(in) :: stock_index
529  integer :: USER_tracer_stock
530 ! This function calculates the mass-weighted integral of all tracer stocks,
531 ! returning the number of stocks it has calculated. If the stock_index
532 ! is present, only the stock corresponding to that coded index is returned.
533 
534 ! Arguments: h - Layer thickness, in m or kg m-2.
535 ! (out) stocks - the mass-weighted integrated amount of each tracer,
536 ! in kg times concentration units.
537 ! (in) G - The ocean's grid structure.
538 ! (in) GV - The ocean's vertical grid structure.
539 ! (in) CS - The control structure returned by a previous call to
540 ! register_USER_tracer.
541 ! (out) names - the names of the stocks calculated.
542 ! (out) units - the units of the stocks calculated.
543 ! (in,opt) stock_index - the coded index of a specific stock being sought.
544 ! Return value: the number of stocks calculated here.
545 
546  integer :: i, j, k, is, ie, js, je, nz, m
547  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
548 
549  user_tracer_stock = 0
550  if (.not.associated(cs)) return
551 
552  if (present(stock_index)) then ; if (stock_index > 0) then
553  ! Check whether this stock is available from this routine.
554 
555  ! No stocks from this routine are being checked yet. Return 0.
556  return
557  endif ; endif
558 
559  do m=1,ntr
560  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="USER_tracer_stock")
561  units(m) = trim(units(m))//" kg"
562  stocks(m) = 0.0
563  do k=1,nz ; do j=js,je ; do i=is,ie
564  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
565  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
566  enddo ; enddo ; enddo
567  stocks(m) = gv%H_to_kg_m2 * stocks(m)
568  enddo
569  user_tracer_stock = ntr
570 
571 end function user_tracer_stock
572 
573 subroutine user_tracer_surface_state(state, h, G, CS)
574  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
575  type(surface), intent(inout) :: state
576  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
577  type(user_tracer_example_cs), pointer :: CS
578 ! This particular tracer package does not report anything back to the coupler.
579 ! The code that is here is just a rough guide for packages that would.
580 ! Arguments: state - A structure containing fields that describe the
581 ! surface state of the ocean.
582 ! (in) h - Layer thickness, in m or kg m-2.
583 ! (in) G - The ocean's grid structure.
584 ! (in) CS - The control structure returned by a previous call to
585 ! register_USER_tracer.
586  integer :: i, j, m, is, ie, js, je
587  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
588 
589  if (.not.associated(cs)) return
590 
591  if (cs%coupled_tracers) then
592  do m=1,ntr
593  ! This call loads the surface vlues into the appropriate array in the
594  ! coupler-type structure.
595  call set_coupler_values(cs%tr(:,:,1,1), state%tr_fields, cs%ind_tr(m), &
596  ind_csurf, is, ie, js, je)
597  enddo
598  endif
599 
600 end subroutine user_tracer_surface_state
601 
602 subroutine user_tracer_example_end(CS)
603  type(user_tracer_example_cs), pointer :: CS
604  integer :: m
605 
606  if (associated(cs)) then
607  if (associated(cs%tr)) deallocate(cs%tr)
608  if (associated(cs%tr_aux)) deallocate(cs%tr_aux)
609  do m=1,ntr
610  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
611  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
612  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
613  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
614  enddo
615 
616  deallocate(cs)
617  endif
618 end subroutine user_tracer_example_end
619 
620 end module user_tracer_example
The following structure contains pointers to various fields which may be used describe the surface st...
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)
subroutine, public user_initialize_tracer(restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp)
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
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.
integer function, public user_tracer_stock(h, stocks, G, GV, CS, names, units, stock_index)
Container for horizontal index ranges for data, computational and global domains. ...
subroutine, public user_tracer_surface_state(state, h, G, CS)
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.
subroutine, public tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS)
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
subroutine, public user_tracer_example_end(CS)
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.
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)
logical function, public user_register_tracer_example(HI, GV, param_file, CS, tr_Reg, restart_CS)
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
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
Controls where open boundary conditions are applied.
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)
integer, parameter ntr