MOM6
RGC_tracer.F90
1 !> This module contains the routines used to set up a
2 !! dynamically passive tracer.
3 !! Set up and use passive tracers requires the following:
4 !! (1) register_RGC_tracer
5 !! (2) apply diffusion, physics/chemistry and advect the tracer
6 
7 !********+*********+*********+*********+*********+*********+*********+**
8 !* *
9 !* By Elizabeth Yankovsky, June 2019 *
10 !*********+*********+*********+*********+*********+*********+***********
11 
12 module rgc_tracer
13 
14 ! This file is part of MOM6. See LICENSE.md for the license.
15 
16 use mom_diag_mediator, only : diag_ctrl
17 use mom_error_handler, only : mom_error, fatal, warning
19 use mom_forcing_type, only : forcing
20 use mom_hor_index, only : hor_index_type
21 use mom_grid, only : ocean_grid_type
22 use mom_io, only : file_exists, read_data, slasher, vardesc, var_desc, query_vardesc
23 use mom_restart, only : mom_restart_cs
24 use mom_ale_sponge, only : set_up_ale_sponge_field, ale_sponge_cs, get_ale_sponge_nz_data
25 use mom_sponge, only : set_up_sponge_field, sponge_cs
26 use mom_time_manager, only : time_type
27 use mom_tracer_registry, only : register_tracer, tracer_registry_type
28 use mom_tracer_diabatic, only : tracer_vertdiff, applytracerboundaryfluxesinout
29 use mom_variables, only : surface
32 
33 implicit none ; private
34 
35 #include <MOM_memory.h>
36 
37 !< Publicly available functions
38 public register_rgc_tracer, initialize_rgc_tracer
39 public rgc_tracer_column_physics, rgc_tracer_end
40 
41 integer, parameter :: ntr = 1 !< The number of tracers in this module.
42 
43 !> tracer control structure
44 type, public :: rgc_tracer_cs ; private
45  logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
46  character(len = 200) :: tracer_ic_file !< The full path to the IC file, or " " to initialize internally.
47  type(time_type), pointer :: time !< A pointer to the ocean model's clock.
48  type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the tracer registry.
49  real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this package.
50  real, pointer :: tr_aux(:,:,:,:) => null() !< The masked tracer concentration.
51  real :: land_val(ntr) = -1.0 !< The value of tr used where land is masked out.
52  real :: lenlat !< the latitudinal or y-direction length of the domain.
53  real :: lenlon !< the longitudinal or x-direction length of the domain.
54  real :: csl !< The length of the continental shelf (x dir, km)
55  real :: lensponge !< the length of the sponge layer.
56  logical :: mask_tracers !< If true, tracers are masked out in massless layers.
57  logical :: use_sponge !< If true, sponges may be applied somewhere in the domain.
58  type(diag_ctrl), pointer :: diag !< A structure that is used to regulate the timing of diagnostic output.
59  type(vardesc) :: tr_desc(ntr) !< Descriptions and metadata for the tracers.
60 end type rgc_tracer_cs
61 
62 contains
63 
64 
65 !> This subroutine is used to register tracer fields
66 function register_rgc_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
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.
149 end function register_rgc_tracer
150 
151 !> Initializes the NTR tracer fields in tr(:,:,:,:)
152 !! and it sets up the tracer output.
153 subroutine initialize_rgc_tracer(restart, day, G, GV, h, diag, OBC, CS, &
154  layer_CSp, sponge_CSp)
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 
273 end subroutine initialize_rgc_tracer
274 
275 !> This subroutine applies diapycnal diffusion and any other column
276 !! tracer physics or chemistry to the tracers from this file.
277 !! This is a simple example of a set of advected passive tracers.
278 subroutine rgc_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
279  evap_CFL_limit, minimum_forcing_depth)
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 
341 end subroutine rgc_tracer_column_physics
342 
343 subroutine rgc_tracer_end(CS)
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
351 end subroutine rgc_tracer_end
352 
353 end module rgc_tracer
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
rgc_tracer
This module contains the routines used to set up a dynamically passive tracer. Set up and use passive...
Definition: RGC_tracer.F90:12
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_ale_sponge::ale_sponge_cs
ALE sponge control structure.
Definition: MOM_ALE_sponge.F90:85
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_hor_index
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Definition: MOM_hor_index.F90:2
mom_ale_sponge
This module contains the routines used to apply sponge layers when using the ALE mode.
Definition: MOM_ALE_sponge.F90:11
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
mom_tracer_diabatic
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
Definition: MOM_tracer_diabatic.F90:4
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_ale_sponge::set_up_ale_sponge_field
Store the reference profile at h points for a variable.
Definition: MOM_ALE_sponge.F90:34
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_sponge
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
mom_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
mom_hor_index::hor_index_type
Container for horizontal index ranges for data, computational and global domains.
Definition: MOM_hor_index.F90:15
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:186
mom_sponge::sponge_cs
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:40
rgc_tracer::rgc_tracer_cs
tracer control structure
Definition: RGC_tracer.F90:44
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:49
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239