MOM6
rgc_tracer Module Reference

Detailed Description

This module contains the routines used to set up a dynamically passive tracer. Set up and use passive tracers requires the following: (1) register_RGC_tracer (2) apply diffusion, physics/chemistry and advect the tracer.

Data Types

type  rgc_tracer_cs
 tracer control structure More...
 

Functions/Subroutines

logical function, public register_rgc_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 This subroutine is used to register tracer fields. More...
 
subroutine, public initialize_rgc_tracer (restart, day, G, GV, h, diag, OBC, CS, layer_CSp, sponge_CSp)
 Initializes the NTR tracer fields in tr(:,:,:,:) and it sets up the tracer output. More...
 
subroutine, public rgc_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, US, 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 rgc_tracer_end (CS)
 

Variables

integer, parameter ntr = 1
 The number of tracers in this module. More...
 

Function/Subroutine Documentation

◆ initialize_rgc_tracer()

subroutine, public rgc_tracer::initialize_rgc_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( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(rgc_tracer_cs), pointer  CS,
type(sponge_cs), pointer  layer_CSp,
type(ale_sponge_cs), pointer  sponge_CSp 
)

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

Parameters
[in]gGrid 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 thickness, in m or kg m-2.
[in]diagStructure used to regulate diagnostic output.
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used. This is not being used for now.
csThe control structure returned by a previous call to RGC_register_tracer.
layer_cspA pointer to the control structure
sponge_cspA pointer to the control structure for the sponges, if they are in use. Otherwise this may be unassociated.

Definition at line 156 of file RGC_tracer.F90.

156 
157  type(ocean_grid_type), intent(in) :: G !< Grid structure.
158  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
159  logical, intent(in) :: restart !< .true. if the fields have already
160  !! been read from a restart file.
161  type(time_type), target, intent(in) :: day !< Time of the start of the run.
162  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
163  intent(in) :: h !< Layer thickness, in m or kg m-2.
164  type(diag_ctrl), target, intent(in) :: diag !< Structure used to regulate diagnostic output.
165  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
166  !! whether, where, and what open boundary
167  !! conditions are used. This is not being used for now.
168  type(RGC_tracer_CS), pointer :: CS !< The control structure returned by a previous
169  !! call to RGC_register_tracer.
170  type(sponge_CS), pointer :: layer_CSp !< A pointer to the control structure
171  type(ALE_sponge_CS), pointer :: sponge_CSp !< A pointer to the control structure for the
172  !! sponges, if they are in use. Otherwise this may be unassociated.
173 
174  real, allocatable :: temp(:,:,:)
175  real, pointer, dimension(:,:,:) :: &
176  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
177  obc_tr1_v => null() ! specify the values of tracer 1 that should come
178  ! in through u- and v- points through the open
179  ! boundary conditions, in the same units as tr.
180  character(len=16) :: name ! A variable's name in a NetCDF file.
181  character(len=72) :: longname ! The long name of that variable.
182  character(len=48) :: units ! The dimensions of the variable.
183  character(len=48) :: flux_units ! The units for tracer fluxes, usually
184  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
185  real, pointer :: tr_ptr(:,:,:) => null()
186  real :: h_neglect ! A thickness that is so small it is usually lost
187  ! in roundoff and can be neglected [H ~> m or kg-2].
188  real :: e(SZK_(G)+1), e_top, e_bot, d_tr ! Heights [Z ~> m].
189  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
190  integer :: IsdB, IedB, JsdB, JedB
191  integer :: nzdata
192 
193  if (.not.associated(cs)) return
194  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
195  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
196  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
197  h_neglect = gv%H_subroundoff
198 
199  cs%Time => day
200  cs%diag => diag
201 
202  if (.not.restart) then
203  if (len_trim(cs%tracer_IC_file) >= 1) then
204  ! Read the tracer concentrations from a netcdf file.
205  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
206  call mom_error(fatal, "RGC_initialize_tracer: Unable to open "// &
207  cs%tracer_IC_file)
208  do m=1,ntr
209  call query_vardesc(cs%tr_desc(m), name, caller="initialize_RGC_tracer")
210  call read_data(cs%tracer_IC_file, trim(name), &
211  cs%tr(:,:,:,m), domain=g%Domain%mpp_domain)
212  enddo
213  else
214  do m=1,ntr
215  do k=1,nz ; do j=js,je ; do i=is,ie
216  cs%tr(i,j,k,m) = 0.0
217  enddo ; enddo ; enddo
218  enddo
219  m=1
220  do j=js,je ; do i=is,ie
221  !set tracer to 1.0 in the surface of the continental shelf
222  if (g%geoLonT(i,j) <= (cs%CSL)) then
223  cs%tr(i,j,1,m) = 1.0 !first layer
224  endif
225  enddo ; enddo
226 
227  endif
228  endif ! restart
229 
230  if ( cs%use_sponge ) then
231 ! If sponges are used, this damps values to zero in the offshore boundary.
232 ! For any tracers that are not damped in the sponge, the call
233 ! to set_up_sponge_field can simply be omitted.
234  if (associated(sponge_csp)) then !ALE mode
235  nzdata = get_ale_sponge_nz_data(sponge_csp)
236  if (nzdata>0) then
237  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nzdata))
238  do k=1,nzdata ; do j=js,je ; do i=is,ie
239  if (g%geoLonT(i,j) >= (cs%lenlon - cs%lensponge) .AND. g%geoLonT(i,j) <= cs%lenlon) then
240  temp(i,j,k) = 0.0
241  endif
242  enddo ; enddo; enddo
243  do m=1,1
244  ! This is needed to force the compiler not to do a copy in the sponge calls.
245  tr_ptr => cs%tr(:,:,:,m)
246  call set_up_ale_sponge_field(temp, g, tr_ptr, sponge_csp)
247  enddo
248  deallocate(temp)
249  endif
250 
251  elseif (associated(layer_csp)) then !layer mode
252  if (nz>0) then
253  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
254  do k=1,nz ; do j=js,je ; do i=is,ie
255  if (g%geoLonT(i,j) >= (cs%lenlon - cs%lensponge) .AND. g%geoLonT(i,j) <= cs%lenlon) then
256  temp(i,j,k) = 0.0
257  endif
258  enddo ; enddo; enddo
259  do m=1,1
260  tr_ptr => cs%tr(:,:,:,m)
261  call set_up_sponge_field(temp, tr_ptr, g, nz, layer_csp)
262  enddo
263  deallocate(temp)
264  endif
265  else
266  call mom_error(fatal, "RGC_initialize_tracer: "// &
267  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
268  endif !selecting mode/calling error if no pointer
269  endif !using sponge
270 

References mom_ale_sponge::get_ale_sponge_nz_data(), mom_error_handler::mom_error(), ntr, mom_io::query_vardesc(), and mom_sponge::set_up_sponge_field().

Here is the call graph for this function:

◆ register_rgc_tracer()

logical function, public rgc_tracer::register_rgc_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(rgc_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

This subroutine is used to register tracer fields.

Parameters
[in]hiA horizontal index type structure.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure indicating the open file to parse for model parameter values.
csA pointer that is set to point to the control structure for this module (in/out).
tr_regA pointer to the tracer registry.
restart_csA pointer to the restart control structure.

Definition at line 68 of file RGC_tracer.F90.

68  type(hor_index_type), intent(in) :: HI !< A horizontal index type structure.
69  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
70  type(param_file_type), intent(in) :: param_file !<A structure indicating the open file to parse
71  !! for model parameter values.
72  type(RGC_tracer_CS), pointer :: CS !< A pointer that is set to point to the control
73  !! structure for this module (in/out).
74  type(tracer_registry_type), pointer :: tr_Reg !< A pointer to the tracer registry.
75  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure.
76 
77  character(len=80) :: name, longname
78 ! This include declares and sets the variable "version".
79 #include "version_variable.h"
80  character(len=40) :: mdl = "RGC_tracer" ! This module's name.
81  character(len=200) :: inputdir
82  real, pointer :: tr_ptr(:,:,:) => null()
83  logical :: register_RGC_tracer
84  integer :: isd, ied, jsd, jed, nz, m
85  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
86 
87  if (associated(cs)) then
88  call mom_error(warning, "RGC_register_tracer called with an "// &
89  "associated control structure.")
90  return
91  endif
92  allocate(cs)
93 
94  ! Read all relevant parameters and write them to the model log.
95  call log_version(param_file, mdl, version, "")
96  call get_param(param_file, mdl, "RGC_TRACER_IC_FILE", cs%tracer_IC_file, &
97  "The name of a file from which to read the initial \n"//&
98  "conditions for the RGC tracers, or blank to initialize \n"//&
99  "them internally.", default=" ")
100  if (len_trim(cs%tracer_IC_file) >= 1) then
101  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
102  inputdir = slasher(inputdir)
103  cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
104  call log_param(param_file, mdl, "INPUTDIR/RGC_TRACER_IC_FILE", &
105  cs%tracer_IC_file)
106  endif
107  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
108  "If true, sponges may be applied anywhere in the domain. \n"//&
109  "The exact location and properties of those sponges are \n"//&
110  "specified from MOM_initialization.F90.", default=.false.)
111 
112  call get_param(param_file, mdl, "LENLAT", cs%lenlat, &
113  "The latitudinal or y-direction length of the domain", &
114  fail_if_missing=.true., do_not_log=.true.)
115 
116  call get_param(param_file, mdl, "LENLON", cs%lenlon, &
117  "The longitudinal or x-direction length of the domain", &
118  fail_if_missing=.true., do_not_log=.true.)
119 
120  call get_param(param_file, mdl, "CONT_SHELF_LENGTH", cs%CSL, &
121  "The length of the continental shelf (x dir, km).", &
122  default=15.0)
123 
124  call get_param(param_file, mdl, "LENSPONGE", cs%lensponge, &
125  "The length of the sponge layer (km).", &
126  default=10.0)
127 
128  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
129  if (cs%mask_tracers) then
130  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
131  endif
132 
133  do m=1,ntr
134  if (m < 10) then ; write(name,'("tr_RGC",I1.1)') m
135  else ; write(name,'("tr_RGC",I2.2)') m ; endif
136  write(longname,'("Concentration of RGC Tracer ",I2.2)') m
137  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
138 
139  ! This is needed to force the compiler not to do a copy in the registration calls.
140  tr_ptr => cs%tr(:,:,:,m)
141  ! Register the tracer for horizontal advection & diffusion.
142  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
143  name=name, longname=longname, units="kg kg-1", &
144  registry_diags=.true., flux_units="kg/s", &
145  restart_cs=restart_cs)
146  enddo
147 
148  cs%tr_Reg => tr_reg
149  register_rgc_tracer = .true.

References mom_error_handler::mom_error(), ntr, mom_tracer_registry::register_tracer(), and mom_io::var_desc().

Here is the call graph for this function:

◆ rgc_tracer_column_physics()

subroutine, public rgc_tracer::rgc_tracer_column_physics ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_old,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_new,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  ea,
real, dimension(szi_(g),szj_(g),szk_(g)), 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(unit_scale_type), intent(in)  US,
type(rgc_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.

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 any possible forcing fields. Unused fields have NULL ptrs.
[in]dtThe amount of time covered by this call [T ~> s].
[in]usA dimensional unit scaling type
csThe control structure returned by a previous call.
[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 [H ~> m or kg m-2].

Definition at line 278 of file RGC_tracer.F90.

278  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
279  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
280  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
281  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
282  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
283  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
284  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
285  intent(in) :: ea !< an array to which the amount of fluid entrained
286  !! from the layer above during this call will be
287  !! added [H ~> m or kg m-2].
288  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
289  intent(in) :: eb !< an array to which the amount of fluid entrained
290  !! from the layer below during this call will be
291  !! added [H ~> m or kg m-2].
292  type(forcing), intent(in) :: fluxes !< A structure containing pointers to any possible
293  !! forcing fields. Unused fields have NULL ptrs.
294  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s].
295  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
296  type(RGC_tracer_CS), pointer :: CS !< The control structure returned by a previous call.
297  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can be
298  !! fluxed out of the top layer in a timestep [nondim].
299  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which fluxes
300  !! can be applied [H ~> m or kg m-2].
301 
302 ! The arguments to this subroutine are redundant in that
303 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
304 
305  real :: b1(SZI_(G)) ! b1 and c1 are variables used by the
306  real :: c1(SZI_(G),SZK_(G)) ! tridiagonal solver.
307  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
308  real :: in_flux(SZI_(G),SZJ_(G),2) ! total amount of tracer to be injected
309 
310  integer :: i, j, k, is, ie, js, je, nz, m
311  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
312 
313  if (.not.associated(cs)) return
314 
315  in_flux(:,:,:) = 0.0
316  m=1
317  do j=js,je ; do i=is,ie
318  !set tracer to 1.0 in the surface of the continental shelf
319  if (g%geoLonT(i,j) <= (cs%CSL)) then
320  cs%tr(i,j,1,m) = 1.0 !first layer
321  endif
322  enddo ; enddo
323 
324  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
325  do m=1,ntr
326  do k=1,nz ;do j=js,je ; do i=is,ie
327  h_work(i,j,k) = h_old(i,j,k)
328  enddo ; enddo ; enddo;
329  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
330  evap_cfl_limit, minimum_forcing_depth, in_flux(:,:,m))
331 
332  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
333  enddo
334  else
335  do m=1,ntr
336  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
337  enddo
338  endif
339 

References mom_tracer_diabatic::applytracerboundaryfluxesinout(), ntr, and mom_tracer_diabatic::tracer_vertdiff().

Here is the call graph for this function:

◆ rgc_tracer_end()

subroutine, public rgc_tracer::rgc_tracer_end ( type(rgc_tracer_cs), pointer  CS)
Parameters
csThe control structure returned by a previous call to RGC_register_tracer.

Definition at line 343 of file RGC_tracer.F90.

343  type(RGC_tracer_CS), pointer :: CS !< The control structure returned by a previous call to RGC_register_tracer.
344  integer :: m
345 
346  if (associated(cs)) then
347  if (associated(cs%tr)) deallocate(cs%tr)
348  deallocate(cs)
349  endif

Variable Documentation

◆ ntr

integer, parameter rgc_tracer::ntr = 1
private

The number of tracers in this module.

Definition at line 42 of file RGC_tracer.F90.

42 integer, parameter :: NTR = 1 !< The number of tracers in this module.

Referenced by initialize_rgc_tracer(), register_rgc_tracer(), and rgc_tracer_column_physics().