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, 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.
 

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 155 of file RGC_tracer.F90.

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

◆ 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 67 of file RGC_tracer.F90.

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

◆ 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(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 [s].
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 [m].

Definition at line 280 of file RGC_tracer.F90.

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

◆ 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 344 of file RGC_tracer.F90.

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