MOM6
DOME_tracer.F90
Go to the documentation of this file.
1 module dome_tracer
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_hor_index, only : hor_index_type
63 use mom_grid, only : ocean_grid_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 
84 
85 ! ntr is the number of tracers in this module.
86 integer, parameter :: ntr = 11
87 
88 type p3d
89  real, dimension(:,:,:), pointer :: p => null()
90 end type p3d
91 
92 type, public :: dome_tracer_cs ; private
93  logical :: coupled_tracers = .false. ! These tracers are not offered to the
94  ! coupler.
95  character(len=200) :: tracer_ic_file ! The full path to the IC file, or " "
96  ! to initialize internally.
97  type(time_type), pointer :: time ! A pointer to the ocean model's clock.
98  type(tracer_registry_type), pointer :: tr_reg => null()
99  real, pointer :: tr(:,:,:,:) => null() ! The array of tracers used in this
100  ! subroutine, in g m-3?
101  real, pointer :: tr_aux(:,:,:,:) => null() ! The masked tracer concentration
102  ! for output, in g m-3.
103  type(p3d), dimension(NTR) :: &
104  tr_adx, &! Tracer zonal advective fluxes in g m-3 m3 s-1.
105  tr_ady, &! Tracer meridional advective fluxes in g m-3 m3 s-1.
106  tr_dfx, &! Tracer zonal diffusive fluxes in g m-3 m3 s-1.
107  tr_dfy ! Tracer meridional diffusive fluxes in g m-3 m3 s-1.
108  real :: land_val(ntr) = -1.0 ! The value of tr used where land is masked out.
109  logical :: mask_tracers ! If true, tracers are masked out in massless layers.
110  logical :: use_sponge
111 
112  integer, dimension(NTR) :: ind_tr ! Indices returned by aof_set_coupler_flux
113  ! if it is used and the surface tracer concentrations are to be
114  ! provided to the coupler.
115 
116  type(diag_ctrl), pointer :: diag ! A structure that is used to regulate the
117  ! timing of diagnostic output.
118  integer, dimension(NTR) :: id_tracer = -1, id_tr_adx = -1, id_tr_ady = -1
119  integer, dimension(NTR) :: id_tr_dfx = -1, id_tr_dfy = -1
120 
121  type(vardesc) :: tr_desc(ntr)
122 end type dome_tracer_cs
123 
124 contains
125 
126 function register_dome_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
127  type(hor_index_type), intent(in) :: HI
128  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
129  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
130  type(dome_tracer_cs), pointer :: CS
131  type(tracer_registry_type), pointer :: tr_Reg
132  type(mom_restart_cs), pointer :: restart_CS
133 ! This subroutine is used to register tracer fields and subroutines
134 ! to be used with MOM.
135 ! Arguments: HI - A horizontal index type structure.
136 ! (in) param_file - A structure indicating the open file to parse for
137 ! model parameter values.
138 ! (in/out) CS - A pointer that is set to point to the control structure
139 ! for this module
140 ! (in/out) tr_Reg - A pointer to the tracer registry.
141 ! (in) restart_CS - A pointer to the restart control structure.
142  character(len=80) :: name, longname
143 ! This include declares and sets the variable "version".
144 #include "version_variable.h"
145  character(len=40) :: mdl = "DOME_tracer" ! This module's name.
146  character(len=200) :: inputdir
147  real, pointer :: tr_ptr(:,:,:) => null()
148  logical :: register_DOME_tracer
149  integer :: isd, ied, jsd, jed, nz, m
150  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
151 
152  if (associated(cs)) then
153  call mom_error(warning, "DOME_register_tracer called with an "// &
154  "associated control structure.")
155  return
156  endif
157  allocate(cs)
158 
159  ! Read all relevant parameters and write them to the model log.
160  call log_version(param_file, mdl, version, "")
161  call get_param(param_file, mdl, "DOME_TRACER_IC_FILE", cs%tracer_IC_file, &
162  "The name of a file from which to read the initial \n"//&
163  "conditions for the DOME tracers, or blank to initialize \n"//&
164  "them internally.", default=" ")
165  if (len_trim(cs%tracer_IC_file) >= 1) then
166  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
167  inputdir = slasher(inputdir)
168  cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
169  call log_param(param_file, mdl, "INPUTDIR/DOME_TRACER_IC_FILE", &
170  cs%tracer_IC_file)
171  endif
172  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
173  "If true, sponges may be applied anywhere in the domain. \n"//&
174  "The exact location and properties of those sponges are \n"//&
175  "specified from MOM_initialization.F90.", default=.false.)
176 
177  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
178  if (cs%mask_tracers) then
179  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
180  endif
181 
182  do m=1,ntr
183  if (m < 10) then ; write(name,'("tr_D",I1.1)') m
184  else ; write(name,'("tr_D",I2.2)') m ; endif
185  write(longname,'("Concentration of DOME Tracer ",I2.2)') m
186  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
187 
188  ! This is needed to force the compiler not to do a copy in the registration
189  ! calls. Curses on the designers and implementers of Fortran90.
190  tr_ptr => cs%tr(:,:,:,m)
191  ! Register the tracer for the restart file.
192  call register_restart_field(tr_ptr, cs%tr_desc(m), .true., restart_cs)
193  ! Register the tracer for horizontal advection & diffusion.
194  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
195  tr_desc_ptr=cs%tr_desc(m))
196 
197  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
198  ! values to the coupler (if any). This is meta-code and its arguments will
199  ! currently (deliberately) give fatal errors if it is used.
200  if (cs%coupled_tracers) &
201  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
202  flux_type=' ', implementation=' ', caller="register_DOME_tracer")
203  enddo
204 
205  cs%tr_Reg => tr_reg
206  register_dome_tracer = .true.
207 end function register_dome_tracer
208 
209 subroutine initialize_dome_tracer(restart, day, G, GV, h, diag, OBC, CS, &
210  sponge_CSp, diag_to_Z_CSp)
211  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
212  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
213  logical, intent(in) :: restart
214  type(time_type), target, intent(in) :: day
215  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
216  type(diag_ctrl), target, intent(in) :: diag
217  type(ocean_obc_type), pointer :: OBC
218  type(dome_tracer_cs), pointer :: CS
219  type(sponge_cs), pointer :: sponge_CSp
220  type(diag_to_z_cs), pointer :: diag_to_Z_CSp
221 ! This subroutine initializes the NTR tracer fields in tr(:,:,:,:)
222 ! and it sets up the tracer output.
223 
224 ! Arguments: restart - .true. if the fields have already been read from
225 ! a restart file.
226 ! (in) day - Time of the start of the run.
227 ! (in) G - The ocean's grid structure.
228 ! (in) GV - The ocean's vertical grid structure.
229 ! (in) h - Layer thickness, in m or kg m-2.
230 ! (in) diag - A structure that is used to regulate diagnostic output.
231 ! (in) OBC - This open boundary condition type specifies whether, where,
232 ! and what open boundary conditions are used.
233 ! (in/out) CS - The control structure returned by a previous call to
234 ! DOME_register_tracer.
235 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
236 ! they are in use. Otherwise this may be unassociated.
237 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
238 ! in depth space.
239  real, allocatable :: temp(:,:,:)
240  real, pointer, dimension(:,:,:) :: &
241  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
242  obc_tr1_v => null() ! specify the values of tracer 1 that should come
243  ! in through u- and v- points through the open
244  ! boundary conditions, in the same units as tr.
245  character(len=16) :: name ! A variable's name in a NetCDF file.
246  character(len=72) :: longname ! The long name of that variable.
247  character(len=48) :: units ! The dimensions of the variable.
248  character(len=48) :: flux_units ! The units for tracer fluxes, usually
249  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
250  real, pointer :: tr_ptr(:,:,:) => null()
251  real :: PI ! 3.1415926... calculated as 4*atan(1)
252  real :: tr_y ! Initial zonally uniform tracer concentrations.
253  real :: dist2 ! The distance squared from a line, in m2.
254  real :: h_neglect ! A thickness that is so small it is usually lost
255  ! in roundoff and can be neglected, in m.
256  real :: e(szk_(g)+1), e_top, e_bot, d_tr
257  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
258  integer :: IsdB, IedB, JsdB, JedB
259 
260  if (.not.associated(cs)) return
261  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
262  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
263  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
264  h_neglect = gv%H_subroundoff
265 
266  cs%Time => day
267  cs%diag => diag
268 
269  if (.not.restart) then
270  if (len_trim(cs%tracer_IC_file) >= 1) then
271  ! Read the tracer concentrations from a netcdf file.
272  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
273  call mom_error(fatal, "DOME_initialize_tracer: Unable to open "// &
274  cs%tracer_IC_file)
275  do m=1,ntr
276  call query_vardesc(cs%tr_desc(m), name, caller="initialize_DOME_tracer")
277  call read_data(cs%tracer_IC_file, trim(name), &
278  cs%tr(:,:,:,m), domain=g%Domain%mpp_domain)
279  enddo
280  else
281  do m=1,ntr
282  do k=1,nz ; do j=js,je ; do i=is,ie
283  cs%tr(i,j,k,m) = 1.0e-20 ! This could just as well be 0.
284  enddo ; enddo ; enddo
285  enddo
286 
287 ! This sets a stripe of tracer across the basin.
288  do m=2,ntr ; do j=js,je ; do i=is,ie
289  tr_y = 0.0
290  if ((m <= 6) .and. (g%geoLatT(i,j) > (300.0+50.0*real(m-1))) .and. &
291  (g%geolatt(i,j) < (350.0+50.0*real(m-1)))) tr_y = 1.0
292  do k=1,nz
293 ! This adds the stripes of tracer to every layer.
294  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + tr_y
295  enddo
296  enddo; enddo; enddo
297 
298  if (ntr > 7) then
299  do j=js,je ; do i=is,ie
300  e(nz+1) = -g%bathyT(i,j)
301  do k=nz,1,-1
302  e(k) = e(k+1) + h(i,j,k)*gv%H_to_m
303  do m=7,ntr
304  e_top = -600.0*real(m-1) + 3000.0
305  e_bot = -600.0*real(m-1) + 2700.0
306  if (e_top < e(k)) then
307  if (e_top < e(k+1)) then ; d_tr = 0.0
308  elseif (e_bot < e(k+1)) then
309  d_tr = (e_top-e(k+1)) / ((h(i,j,k)+h_neglect)*gv%H_to_m)
310  else ; d_tr = (e_top-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_m)
311  endif
312  elseif (e_bot < e(k)) then
313  if (e_bot < e(k+1)) then ; d_tr = 1.0
314  else ; d_tr = (e(k)-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_m)
315  endif
316  else
317  d_tr = 0.0
318  endif
319  if (h(i,j,k) < 2.0*gv%Angstrom) d_tr=0.0
320  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + d_tr
321  enddo
322  enddo
323  enddo ; enddo
324  endif
325 
326  endif
327  endif ! restart
328 
329  if ( cs%use_sponge ) then
330 ! If sponges are used, this example damps tracers in sponges in the
331 ! northern half of the domain to 1 and tracers in the southern half
332 ! to 0. For any tracers that are not damped in the sponge, the call
333 ! to set_up_sponge_field can simply be omitted.
334  if (.not.associated(sponge_csp)) &
335  call mom_error(fatal, "DOME_initialize_tracer: "// &
336  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
337 
338  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
339  do k=1,nz ; do j=js,je ; do i=is,ie
340  if (g%geoLatT(i,j) > 700.0 .and. (k > nz/2)) then
341  temp(i,j,k) = 1.0
342  else
343  temp(i,j,k) = 0.0
344  endif
345  enddo ; enddo ; enddo
346 
347 ! do m=1,NTR
348  do m=1,1
349  ! This is needed to force the compiler not to do a copy in the sponge
350  ! calls. Curses on the designers and implementers of Fortran90.
351  tr_ptr => cs%tr(:,:,:,m)
352  call set_up_sponge_field(temp, tr_ptr, g, nz, sponge_csp)
353  enddo
354  deallocate(temp)
355  endif
356 
357  if (associated(obc)) then
358  call query_vardesc(cs%tr_desc(1), name, caller="initialize_DOME_tracer")
359  if (obc%specified_v_BCs_exist_globally) then
360  allocate(obc_tr1_v(g%isd:g%ied,g%jsd:g%jed,nz))
361  do k=1,nz ; do j=g%jsd,g%jed ; do i=g%isd,g%ied
362  if (k < nz/2) then ; obc_tr1_v(i,j,k) = 0.0
363  else ; obc_tr1_v(i,j,k) = 1.0 ; endif
364  enddo ; enddo ; enddo
365  call add_tracer_obc_values(trim(name), cs%tr_Reg, &
366  0.0, obc_in_v=obc_tr1_v)
367  else
368  ! This is not expected in the DOME example.
369  call add_tracer_obc_values(trim(name), cs%tr_Reg, 0.0)
370  endif
371  ! All tracers but the first have 0 concentration in their inflows. As this
372  ! is the default value, the following calls are unnecessary.
373  do m=2,ntr
374  call query_vardesc(cs%tr_desc(m), name, caller="initialize_DOME_tracer")
375  call add_tracer_obc_values(trim(name), cs%tr_Reg, 0.0)
376  enddo
377  endif
378 
379  ! This needs to be changed if the units of tracer are changed above.
380  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
381  else ; flux_units = "kg s-1" ; endif
382 
383  do m=1,ntr
384  ! Register the tracer for the restart file.
385  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
386  caller="initialize_DOME_tracer")
387  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
388  day, trim(longname) , trim(units))
389  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
390  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
391  trim(flux_units))
392  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
393  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
394  trim(flux_units))
395  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
396  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
397  trim(flux_units))
398  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
399  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
400  trim(flux_units))
401  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
402  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
403  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
404  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
405 
406 ! Register the tracer for horizontal advection & diffusion.
407  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
408  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
409  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
410  cs%tr_ady(m)%p, cs%tr_dfx(m)%p, cs%tr_dfy(m)%p)
411 
412  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
413  day, g, diag_to_z_csp)
414  enddo
415 
416 end subroutine initialize_dome_tracer
417 
418 subroutine dome_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
419  evap_CFL_limit, minimum_forcing_depth)
420  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
421  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
422  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
423  type(forcing), intent(in) :: fluxes
424  real, intent(in) :: dt !< The amount of time covered by this call, in s
425  type(dome_tracer_cs), pointer :: CS
426  real, optional,intent(in) :: evap_CFL_limit
427  real, optional,intent(in) :: minimum_forcing_depth
428 ! This subroutine applies diapycnal diffusion and any other column
429 ! tracer physics or chemistry to the tracers from this file.
430 ! This is a simple example of a set of advected passive tracers.
431 
432 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
433 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
434 ! (in) ea - an array to which the amount of fluid entrained
435 ! from the layer above during this call will be
436 ! added, in m or kg m-2.
437 ! (in) eb - an array to which the amount of fluid entrained
438 ! from the layer below during this call will be
439 ! added, in m or kg m-2.
440 ! (in) fluxes - A structure containing pointers to any possible
441 ! forcing fields. Unused fields have NULL ptrs.
442 ! (in) dt - The amount of time covered by this call, in s.
443 ! (in) G - The ocean's grid structure.
444 ! (in) GV - The ocean's vertical grid structure.
445 ! (in) CS - The control structure returned by a previous call to
446 ! DOME_register_tracer.
447 !
448 ! The arguments to this subroutine are redundant in that
449 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
450 
451  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
452  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
453  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
454  integer :: i, j, k, is, ie, js, je, nz, m
455  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
456 
457  if (.not.associated(cs)) return
458 
459  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
460  do m=1,ntr
461  do k=1,nz ;do j=js,je ; do i=is,ie
462  h_work(i,j,k) = h_old(i,j,k)
463  enddo ; enddo ; enddo;
464  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
465  evap_cfl_limit, minimum_forcing_depth)
466  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
467  enddo
468  else
469  do m=1,ntr
470  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
471  enddo
472  endif
473 
474  if (cs%mask_tracers) then
475  do m = 1,ntr ; if (cs%id_tracer(m) > 0) then
476  do k=1,nz ; do j=js,je ; do i=is,ie
477  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
478  cs%tr_aux(i,j,k,m) = cs%land_val(m)
479  else
480  cs%tr_aux(i,j,k,m) = cs%tr(i,j,k,m)
481  endif
482  enddo ; enddo ; enddo
483  endif ; enddo
484  endif
485 
486  do m=1,ntr
487  if (cs%mask_tracers) then
488  if (cs%id_tracer(m)>0) &
489  call post_data(cs%id_tracer(m),cs%tr_aux(:,:,:,m),cs%diag)
490  else
491  if (cs%id_tracer(m)>0) &
492  call post_data(cs%id_tracer(m),cs%tr(:,:,:,m),cs%diag)
493  endif
494  if (cs%id_tr_adx(m)>0) &
495  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
496  if (cs%id_tr_ady(m)>0) &
497  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
498  if (cs%id_tr_dfx(m)>0) &
499  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
500  if (cs%id_tr_dfy(m)>0) &
501  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
502  enddo
503 
504 end subroutine dome_tracer_column_physics
505 
506 subroutine dome_tracer_surface_state(state, h, G, CS)
507  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
508  type(surface), intent(inout) :: state
509  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
510  type(dome_tracer_cs), pointer :: CS
511 ! This particular tracer package does not report anything back to the coupler.
512 ! The code that is here is just a rough guide for packages that would.
513 ! Arguments: state - A structure containing fields that describe the
514 ! surface state of the ocean.
515 ! (in) h - Layer thickness, in m or kg m-2.
516 ! (in) G - The ocean's grid structure.
517 ! (in) CS - The control structure returned by a previous call to
518 ! DOME_register_tracer.
519  integer :: i, j, m, is, ie, js, je
520  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
521 
522  if (.not.associated(cs)) return
523 
524  if (cs%coupled_tracers) then
525  do m=1,ntr
526  ! This call loads the surface vlues into the appropriate array in the
527  ! coupler-type structure.
528  call set_coupler_values(cs%tr(:,:,1,1), state%tr_fields, cs%ind_tr(m), &
529  ind_csurf, is, ie, js, je)
530  enddo
531  endif
532 
533 end subroutine dome_tracer_surface_state
534 
535 subroutine dome_tracer_end(CS)
536  type(dome_tracer_cs), pointer :: CS
537  integer :: m
538 
539  if (associated(cs)) then
540  if (associated(cs%tr)) deallocate(cs%tr)
541  if (associated(cs%tr_aux)) deallocate(cs%tr_aux)
542  do m=1,ntr
543  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
544  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
545  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
546  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
547  enddo
548 
549  deallocate(cs)
550  endif
551 end subroutine dome_tracer_end
552 
553 end module dome_tracer
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)
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
subroutine, public dome_tracer_surface_state(state, h, G, CS)
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 dome_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
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...
subroutine, public initialize_dome_tracer(restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp)
logical function, public register_dome_tracer(HI, GV, param_file, CS, tr_Reg, restart_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
Controls where open boundary conditions are applied.
subroutine, public dome_tracer_end(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)
integer, parameter ntr
Definition: DOME_tracer.F90:86