MOM6
dyed_obc_tracer Module Reference

Detailed Description

This tracer package dyes flow through open boundaries.

By Kate Hedstrom, 2017, copied from DOME tracers and also dye_example.

This file contains an example of the code that is needed to set up and use a set of dynamically passive tracers. These tracers dye the inflowing water, one per open boundary segment.

A single subroutine is called from within each file to register each of the tracers for reinitialization and advection and to register the subroutine that initializes the tracers and set up their output and the subroutine that does any tracer physics or chemistry along with diapycnal mixing (included here because some tracers may float or swim vertically or dye diapycnal processes).

Data Types

type  dyed_obc_tracer_cs
 The control structure for the dyed_obc tracer package. More...
 

Functions/Subroutines

logical function, public register_dyed_obc_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Register tracer fields and subroutines to be used with MOM. More...
 
subroutine, public initialize_dyed_obc_tracer (restart, day, G, GV, h, diag, OBC, CS)
 Initializes the CSntr tracer fields in tr(:,:,:,:) and sets up the tracer output. More...
 
subroutine, public dyed_obc_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
 This subroutine applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this file. This is a simple example of a set of advected passive tracers. More...
 
subroutine, public dyed_obc_tracer_end (CS)
 Clean up memory allocations, if any. More...
 

Function/Subroutine Documentation

◆ dyed_obc_tracer_column_physics()

subroutine, public dyed_obc_tracer::dyed_obc_tracer_column_physics ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_old,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_new,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  ea,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  eb,
type(forcing), intent(in)  fluxes,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(dyed_obc_tracer_cs), pointer  CS,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

This subroutine applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this file. This is a simple example of a set of advected passive tracers.

The arguments to this subroutine are redundant in that h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]h_oldLayer thickness before entrainment [H ~> m or kg m-2].
[in]h_newLayer thickness after entrainment [H ~> m or kg m-2].
[in]eaan array to which the amount of fluid entrained
[in]eban array to which the amount of fluid entrained
[in]fluxesA structure containing pointers to thermodynamic and tracer forcing fields. Unused fields have NULL ptrs.
[in]dtThe amount of time covered by this call [s]
csThe control structure returned by a previous call to dyed_obc_register_tracer.
[in]evap_cfl_limitLimit on the fraction of the water that can be fluxed out of the top layer in a timestep [nondim]
[in]minimum_forcing_depthThe smallest depth over which fluxes can be applied [m]

Definition at line 204 of file dyed_obc_tracer.F90.

204  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
205  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
206  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
207  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
208  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
209  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
210  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
211  intent(in) :: ea !< an array to which the amount of fluid entrained
212  !! from the layer above during this call will be
213  !! added [H ~> m or kg m-2].
214  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
215  intent(in) :: eb !< an array to which the amount of fluid entrained
216  !! from the layer below during this call will be
217  !! added [H ~> m or kg m-2].
218  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
219  !! and tracer forcing fields. Unused fields have NULL ptrs.
220  real, intent(in) :: dt !< The amount of time covered by this call [s]
221  type(dyed_obc_tracer_CS), pointer :: CS !< The control structure returned by a previous
222  !! call to dyed_obc_register_tracer.
223  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
224  !! be fluxed out of the top layer in a timestep [nondim]
225  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
226  !! fluxes can be applied [m]
227 
228 ! Local variables
229  real :: b1(SZI_(G)) ! b1 and c1 are variables used by the
230  real :: c1(SZI_(G),SZK_(G)) ! tridiagonal solver.
231  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
232  integer :: i, j, k, is, ie, js, je, nz, m
233  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
234 
235  if (.not.associated(cs)) return
236  if (cs%ntr < 1) return
237 
238  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
239  do m=1,cs%ntr
240  do k=1,nz ;do j=js,je ; do i=is,ie
241  h_work(i,j,k) = h_old(i,j,k)
242  enddo ; enddo ; enddo
243  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
244  evap_cfl_limit, minimum_forcing_depth)
245  if (nz > 1) call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
246  enddo
247  else
248  do m=1,cs%ntr
249  if (nz > 1) call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
250  enddo
251  endif
252 

◆ dyed_obc_tracer_end()

subroutine, public dyed_obc_tracer::dyed_obc_tracer_end ( type(dyed_obc_tracer_cs), pointer  CS)

Clean up memory allocations, if any.

Parameters
csThe control structure returned by a previous call to dyed_obc_register_tracer.

Definition at line 257 of file dyed_obc_tracer.F90.

257  type(dyed_obc_tracer_CS), pointer :: CS !< The control structure returned by a previous
258  !! call to dyed_obc_register_tracer.
259  integer :: m
260 
261  if (associated(cs)) then
262  if (associated(cs%tr)) deallocate(cs%tr)
263 
264  deallocate(cs)
265  endif

◆ initialize_dyed_obc_tracer()

subroutine, public dyed_obc_tracer::initialize_dyed_obc_tracer ( logical, intent(in)  restart,
type(time_type), intent(in), target  day,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(dyed_obc_tracer_cs), pointer  CS 
)

Initializes the CSntr tracer fields in tr(:,:,:,:) and sets up the tracer output.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]diagStructure used to regulate diagnostic output.
obcStructure specifying open boundary options.
csThe control structure returned by a previous call to dyed_obc_register_tracer.

Definition at line 135 of file dyed_obc_tracer.F90.

135  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
136  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
137  logical, intent(in) :: restart !< .true. if the fields have already
138  !! been read from a restart file.
139  type(time_type), target, intent(in) :: day !< Time of the start of the run.
140  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
141  type(diag_ctrl), target, intent(in) :: diag !< Structure used to regulate diagnostic output.
142  type(ocean_OBC_type), pointer :: OBC !< Structure specifying open boundary options.
143  type(dyed_obc_tracer_CS), pointer :: CS !< The control structure returned by a previous
144  !! call to dyed_obc_register_tracer.
145 
146 ! Local variables
147  real, allocatable :: temp(:,:,:)
148  real, pointer, dimension(:,:,:) :: &
149  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
150  obc_tr1_v => null() ! specify the values of tracer 1 that should come
151  ! in through u- and v- points through the open
152  ! boundary conditions, in the same units as tr.
153  character(len=24) :: name ! A variable's name in a NetCDF file.
154  character(len=72) :: longname ! The long name of that variable.
155  character(len=48) :: units ! The dimensions of the variable.
156  character(len=48) :: flux_units ! The units for tracer fluxes, usually
157  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
158  real, pointer :: tr_ptr(:,:,:) => null()
159  real :: h_neglect ! A thickness that is so small it is usually lost
160  ! in roundoff and can be neglected [H ~> m or kg m-2].
161  real :: e(SZK_(G)+1), e_top, e_bot, d_tr
162  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
163  integer :: IsdB, IedB, JsdB, JedB
164 
165  if (.not.associated(cs)) return
166  if (cs%ntr < 1) return
167  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
168  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
169  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
170  h_neglect = gv%H_subroundoff
171 
172  cs%Time => day
173  cs%diag => diag
174 
175  if (.not.restart) then
176  if (len_trim(cs%tracer_IC_file) >= 1) then
177  ! Read the tracer concentrations from a netcdf file.
178  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
179  call mom_error(fatal, "dyed_obc_initialize_tracer: Unable to open "// &
180  cs%tracer_IC_file)
181  do m=1,cs%ntr
182  call query_vardesc(cs%tr_desc(m), name, caller="initialize_dyed_obc_tracer")
183  call mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
184  enddo
185  else
186  do m=1,cs%ntr
187  do k=1,nz ; do j=js,je ; do i=is,ie
188  cs%tr(i,j,k,m) = 0.0
189  enddo ; enddo ; enddo
190  enddo
191  endif
192  endif ! restart
193 

◆ register_dyed_obc_tracer()

logical function, public dyed_obc_tracer::register_dyed_obc_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(dyed_obc_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

Register tracer fields and subroutines to be used with MOM.

Parameters
[in]hiA horizontal index type structure.
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters
csA pointer that is set to point to the control structure for this module
tr_regA pointer to the tracer registry.
restart_csA pointer to the restart control structure.

Definition at line 54 of file dyed_obc_tracer.F90.

54  type(hor_index_type), intent(in) :: HI !< A horizontal index type structure.
55  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
56  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
57  type(dyed_obc_tracer_CS), pointer :: CS !< A pointer that is set to point to the
58  !! control structure for this module
59  type(tracer_registry_type), pointer :: tr_Reg !< A pointer to the tracer registry.
60  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure.
61 
62 ! Local variables
63  character(len=80) :: name, longname
64 ! This include declares and sets the variable "version".
65 #include "version_variable.h"
66  character(len=40) :: mdl = "dyed_obc_tracer" ! This module's name.
67  character(len=200) :: inputdir
68  character(len=48) :: flux_units ! The units for tracer fluxes, usually
69  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
70  real, pointer :: tr_ptr(:,:,:) => null()
71  logical :: register_dyed_obc_tracer
72  integer :: isd, ied, jsd, jed, nz, m
73  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
74 
75  if (associated(cs)) then
76  call mom_error(warning, "dyed_obc_register_tracer called with an "// &
77  "associated control structure.")
78  return
79  endif
80  allocate(cs)
81 
82  ! Read all relevant parameters and write them to the model log.
83  call log_version(param_file, mdl, version, "")
84  call get_param(param_file, mdl, "NUM_DYE_TRACERS", cs%ntr, &
85  "The number of dye tracers in this run. Each tracer "//&
86  "should have a separate boundary segment.", default=0)
87  allocate(cs%ind_tr(cs%ntr))
88  allocate(cs%tr_desc(cs%ntr))
89 
90  call get_param(param_file, mdl, "dyed_obc_TRACER_IC_FILE", cs%tracer_IC_file, &
91  "The name of a file from which to read the initial "//&
92  "conditions for the dyed_obc tracers, or blank to initialize "//&
93  "them internally.", default=" ")
94  if (len_trim(cs%tracer_IC_file) >= 1) then
95  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
96  inputdir = slasher(inputdir)
97  cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
98  call log_param(param_file, mdl, "INPUTDIR/dyed_obc_TRACER_IC_FILE", &
99  cs%tracer_IC_file)
100  endif
101 
102  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
103 
104  do m=1,cs%ntr
105  write(name,'("dye_",I2.2)') m
106  write(longname,'("Concentration of dyed_obc Tracer ",I2.2)') m
107  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
108  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
109  else ; flux_units = "kg s-1" ; endif
110 
111  ! This is needed to force the compiler not to do a copy in the registration
112  ! calls. Curses on the designers and implementers of Fortran90.
113  tr_ptr => cs%tr(:,:,:,m)
114  ! Register the tracer for horizontal advection, diffusion, and restarts.
115  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
116  name=name, longname=longname, units="kg kg-1", &
117  registry_diags=.true., flux_units=flux_units, &
118  restart_cs=restart_cs)
119 
120  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
121  ! values to the coupler (if any). This is meta-code and its arguments will
122  ! currently (deliberately) give fatal errors if it is used.
123  if (cs%coupled_tracers) &
124  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
125  flux_type=' ', implementation=' ', caller="register_dyed_obc_tracer")
126  enddo
127 
128  cs%tr_Reg => tr_reg
129  cs%restart_CSp => restart_cs
130  register_dyed_obc_tracer = .true.