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, mom_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
30 use mom_variables, only : surface
33 
34 implicit none ; private
35 
36 #include <MOM_memory.h>
37 
38 !< Publicly available functions
39 public register_rgc_tracer, initialize_rgc_tracer
40 public rgc_tracer_column_physics, rgc_tracer_end
41 
42 integer, parameter :: ntr = 1 !< The number of tracers in this module.
43 
44 !> tracer control structure
45 type, public :: rgc_tracer_cs ; private
46  logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
47  character(len = 200) :: tracer_ic_file !< The full path to the IC file, or " " to initialize internally.
48  type(time_type), pointer :: time !< A pointer to the ocean model's clock.
49  type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the tracer registry.
50  real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this package.
51  real, pointer :: tr_aux(:,:,:,:) => null() !< The masked tracer concentration.
52  real :: land_val(ntr) = -1.0 !< The value of tr used where land is masked out.
53  real :: lenlat !< the latitudinal or y-direction length of the domain.
54  real :: lenlon !< the longitudinal or x-direction length of the domain.
55  real :: csl !< The length of the continental shelf (x dir, km)
56  real :: lensponge !< the length of the sponge layer.
57  logical :: mask_tracers !< If true, tracers are masked out in massless layers.
58  logical :: use_sponge !< If true, sponges may be applied somewhere in the domain.
59  type(diag_ctrl), pointer :: diag !< A structure that is used to regulate the timing of diagnostic output.
60  type(vardesc) :: tr_desc(ntr) !< Descriptions and metadata for the tracers.
61 end type rgc_tracer_cs
62 
63 contains
64 
65 
66 !> This subroutine is used to register tracer fields
67 function register_rgc_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
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.
150 end function register_rgc_tracer
151 
152 !> Initializes the NTR tracer fields in tr(:,:,:,:)
153 !! and it sets up the tracer output.
154 subroutine initialize_rgc_tracer(restart, day, G, GV, h, diag, OBC, CS, &
155  layer_CSp, sponge_CSp)
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 mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
211  enddo
212  else
213  do m=1,ntr
214  do k=1,nz ; do j=js,je ; do i=is,ie
215  cs%tr(i,j,k,m) = 0.0
216  enddo ; enddo ; enddo
217  enddo
218  m=1
219  do j=js,je ; do i=is,ie
220  !set tracer to 1.0 in the surface of the continental shelf
221  if (g%geoLonT(i,j) <= (cs%CSL)) then
222  cs%tr(i,j,1,m) = 1.0 !first layer
223  endif
224  enddo ; enddo
225 
226  endif
227  endif ! restart
228 
229  if ( cs%use_sponge ) then
230 ! If sponges are used, this damps values to zero in the offshore boundary.
231 ! For any tracers that are not damped in the sponge, the call
232 ! to set_up_sponge_field can simply be omitted.
233  if (associated(sponge_csp)) then !ALE mode
234  nzdata = get_ale_sponge_nz_data(sponge_csp)
235  if (nzdata>0) then
236  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nzdata))
237  do k=1,nzdata ; do j=js,je ; do i=is,ie
238  if (g%geoLonT(i,j) >= (cs%lenlon - cs%lensponge) .AND. g%geoLonT(i,j) <= cs%lenlon) then
239  temp(i,j,k) = 0.0
240  endif
241  enddo ; enddo; enddo
242  do m=1,1
243  ! This is needed to force the compiler not to do a copy in the sponge calls.
244  tr_ptr => cs%tr(:,:,:,m)
245  call set_up_ale_sponge_field(temp, g, tr_ptr, sponge_csp)
246  enddo
247  deallocate(temp)
248  endif
249 
250  elseif (associated(layer_csp)) then !layer mode
251  if (nz>0) then
252  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
253  do k=1,nz ; do j=js,je ; do i=is,ie
254  if (g%geoLonT(i,j) >= (cs%lenlon - cs%lensponge) .AND. g%geoLonT(i,j) <= cs%lenlon) then
255  temp(i,j,k) = 0.0
256  endif
257  enddo ; enddo; enddo
258  do m=1,1
259  tr_ptr => cs%tr(:,:,:,m)
260  call set_up_sponge_field(temp, tr_ptr, g, nz, layer_csp)
261  enddo
262  deallocate(temp)
263  endif
264  else
265  call mom_error(fatal, "RGC_initialize_tracer: "// &
266  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
267  endif !selecting mode/calling error if no pointer
268  endif !using sponge
269 
270 end subroutine initialize_rgc_tracer
271 
272 !> This subroutine applies diapycnal diffusion and any other column
273 !! tracer physics or chemistry to the tracers from this file.
274 !! This is a simple example of a set of advected passive tracers.
275 subroutine rgc_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, &
276  evap_CFL_limit, minimum_forcing_depth)
277  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
278  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
279  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
280  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
281  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
282  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
283  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
284  intent(in) :: ea !< an array to which the amount of fluid entrained
285  !! from the layer above during this call will be
286  !! added [H ~> m or kg m-2].
287  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
288  intent(in) :: eb !< an array to which the amount of fluid entrained
289  !! from the layer below during this call will be
290  !! added [H ~> m or kg m-2].
291  type(forcing), intent(in) :: fluxes !< A structure containing pointers to any possible
292  !! forcing fields. Unused fields have NULL ptrs.
293  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s].
294  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
295  type(rgc_tracer_cs), pointer :: cs !< The control structure returned by a previous call.
296  real, optional, intent(in) :: evap_cfl_limit !< Limit on the fraction of the water that can be
297  !! fluxed out of the top layer in a timestep [nondim].
298  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which fluxes
299  !! can be applied [H ~> m or kg m-2].
300 
301 ! The arguments to this subroutine are redundant in that
302 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
303 
304  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
305  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
306  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
307  real :: in_flux(szi_(g),szj_(g),2) ! total amount of tracer to be injected
308 
309  integer :: i, j, k, is, ie, js, je, nz, m
310  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
311 
312  if (.not.associated(cs)) return
313 
314  in_flux(:,:,:) = 0.0
315  m=1
316  do j=js,je ; do i=is,ie
317  !set tracer to 1.0 in the surface of the continental shelf
318  if (g%geoLonT(i,j) <= (cs%CSL)) then
319  cs%tr(i,j,1,m) = 1.0 !first layer
320  endif
321  enddo ; enddo
322 
323  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
324  do m=1,ntr
325  do k=1,nz ;do j=js,je ; do i=is,ie
326  h_work(i,j,k) = h_old(i,j,k)
327  enddo ; enddo ; enddo;
328  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
329  evap_cfl_limit, minimum_forcing_depth, in_flux(:,:,m))
330 
331  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
332  enddo
333  else
334  do m=1,ntr
335  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
336  enddo
337  endif
338 
339 end subroutine rgc_tracer_column_physics
340 
341 subroutine rgc_tracer_end(CS)
342  type(rgc_tracer_cs), pointer :: cs !< The control structure returned by a previous call to RGC_register_tracer.
343  integer :: m
344 
345  if (associated(cs)) then
346  if (associated(cs%tr)) deallocate(cs%tr)
347  deallocate(cs)
348  endif
349 end subroutine rgc_tracer_end
350 
351 end module 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
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Wraps the FMS time manager functions.
This module implements boundary forcing for MOM6.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Store the reference profile at h points for a variable.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
Defines the horizontal index type (hor_index_type) used for providing index ranges.
This module contains the routines used to apply sponge layers when using the ALE mode.
An overloaded interface to log the values of various types of parameters.
ALE sponge control structure.
Container for horizontal index ranges for data, computational and global domains. ...
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:75
Describes various unit conversion factors.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Type to carry basic tracer information.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Routines for error handling and I/O management.
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:41
The MOM6 facility for reading and writing restart files, and querying what has been read...
Definition: MOM_restart.F90:2
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
tracer control structure
Definition: RGC_tracer.F90:45
Controls where open boundary conditions are applied.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Read a data field from a file.
Definition: MOM_io.F90:74
An overloaded interface to read and log the values of various types of parameters.