MOM6
MOM_tracer_registry.F90
Go to the documentation of this file.
1 !> This module contains the tracer_registry_type and the subroutines
2 !! that handle registration of tracers and related subroutines.
3 !! The primary subroutine, register_tracer, is called to indicate the
4 !! tracers advected and diffused.
6 
7 ! This file is part of MOM6. See LICENSE.md for the license.
8 
9 ! use MOM_diag_mediator, only : diag_ctrl
10 use mom_coms, only : reproducing_sum
11 use mom_debugging, only : hchksum
12 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg, is_root_pe
14 use mom_hor_index, only : hor_index_type
15 use mom_grid, only : ocean_grid_type
16 use mom_io, only : vardesc, query_vardesc
17 use mom_time_manager, only : time_type
19 
20 implicit none ; private
21 
22 #include <MOM_memory.h>
23 
24 public register_tracer
26 public mom_tracer_chksum
27 public mom_tracer_chkinv
32 
33 !> The tracer type
34 type, public :: tracer_type
35 
36  real, dimension(:,:,:), pointer :: t => null() !< tracer concentration array
37  real :: obc_inflow_conc= 0.0 !< tracer concentration for generic inflows
38  real, dimension(:,:,:), pointer :: obc_in_u => null() !< structured values for flow into the domain
39  !! specified in OBCs through u-face of cell
40  real, dimension(:,:,:), pointer :: obc_in_v => null() !< structured values for flow into the domain
41  !! specified in OBCs through v-face of cell
42 
43  real, dimension(:,:,:), pointer :: ad_x => null() !< diagnostic array for x-advective tracer flux
44  real, dimension(:,:,:), pointer :: ad_y => null() !< diagnostic array for y-advective tracer flux
45  real, dimension(:,:), pointer :: ad2d_x => null() !< diagnostic vertical sum x-advective tracer flux
46  !! in units of (conc * m3/s or conc * kg/s)
47  real, dimension(:,:), pointer :: ad2d_y => null() !< diagnostic vertical sum y-advective tracer flux
48  !! in units of (conc * m3/s or conc * kg/s)
49 
50  real, dimension(:,:,:), pointer :: df_x => null() !< diagnostic array for x-diffusive tracer flux
51  real, dimension(:,:,:), pointer :: df_y => null() !< diagnostic array for y-diffusive tracer flux
52  real, dimension(:,:), pointer :: df2d_x => null() !< diagnostic vertical sum x-diffusive flux
53  !! in units of (conc * m3/s or conc * kg/s)
54  real, dimension(:,:), pointer :: df2d_y => null() !< diagnostic vertical sum y-diffusive flux
55  !! in units of (conc * m3/s or conc * kg/s)
56 
57  real, dimension(:,:,:), pointer :: advection_xy => null() !< convergence of lateral advective tracer fluxes
58 
59  character(len=32) :: name !< tracer name used for error messages
60  type(vardesc), pointer :: vd => null() !< metadata describing the tracer
61 
62 end type tracer_type
63 
64 !> Type to carry basic tracer information
65 type, public :: tracer_registry_type
66  integer :: ntr = 0 !< number of registered tracers
67  type(tracer_type) :: tr(max_fields_) !< array of registered tracers
68 ! type(diag_ctrl), pointer :: diag !< structure to regulate timing of diagnostics
69  logical :: locked = .false. !< New tracers may be registered if locked=.false.
70  !! When locked=.true.,no more tracers can be registered,
71  !! at which point common diagnostics can be set up
72  !! for the registered tracers.
73 end type tracer_registry_type
74 
75 contains
76 
77 !> This subroutine registers a tracer to be advected and laterally diffused.
78 subroutine register_tracer(tr1, tr_desc, param_file, HI, GV, Reg, tr_desc_ptr, ad_x, ad_y,&
79  df_x, df_y, OBC_inflow, OBC_in_u, OBC_in_v, &
80  ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
81  type(hor_index_type), intent(in) :: HI !< horizontal index type
82  type(verticalgrid_type), intent(in) :: GV !< ocean vertical grid structure
83  real, dimension(SZI_(HI),SZJ_(HI),SZK_(GV)), target :: tr1 !< pointer to the tracer (concentration units)
84  type(vardesc), intent(in) :: tr_desc !< metadata about the tracer
85  type(param_file_type), intent(in) :: param_file !< file to parse for model parameter values
86  type(tracer_registry_type), pointer :: Reg !< pointer to the tracer registry
87  type(vardesc), target, optional :: tr_desc_ptr !< A target that can be used to set a pointer to the
88  !! stored value of tr%tr_desc. This target must be an
89  !! enduring part of the control structure, because the tracer
90  !! registry will use this memory, but it also means that any
91  !! updates to this structure in the calling module will be
92  !! available subsequently to the tracer registry.
93  real, pointer, dimension(:,:,:), optional :: ad_x !< diagnostic x-advective flux (CONC m3/s or CONC*kg/s)
94  real, pointer, dimension(:,:,:), optional :: ad_y !< diagnostic y-advective flux (CONC m3/s or CONC*kg/s)
95  real, pointer, dimension(:,:,:), optional :: df_x !< diagnostic x-diffusive flux (CONC m3/s or CONC*kg/s)
96  real, pointer, dimension(:,:,:), optional :: df_y !< diagnostic y-diffusive flux (CONC m3/s or CONC*kg/s)
97 
98  real, intent(in), optional :: OBC_inflow !< the tracer for all inflows via OBC for which OBC_in_u
99  !! or OBC_in_v are not specified (units of tracer CONC)
100  real, pointer, dimension(:,:,:), optional :: OBC_in_u !< tracer at inflows through u-faces of
101  !! tracer cells (units of tracer CONC)
102  real, pointer, dimension(:,:,:), optional :: OBC_in_v !< tracer at inflows through v-faces of
103  !! tracer cells (units of tracer CONC)
104 
105  real, dimension(:,:), pointer, optional :: ad_2d_x !< vert sum of diagnostic x-advect flux (CONC m3/s or CONC*kg/s)
106  real, dimension(:,:), pointer, optional :: ad_2d_y !< vert sum of diagnostic y-advect flux (CONC m3/s or CONC*kg/s)
107  real, dimension(:,:), pointer, optional :: df_2d_x !< vert sum of diagnostic x-diffuse flux (CONC m3/s or CONC*kg/s)
108  real, dimension(:,:), pointer, optional :: df_2d_y !< vert sum of diagnostic y-diffuse flux (CONC m3/s or CONC*kg/s)
109 
110  real, pointer, dimension(:,:,:), optional :: advection_xy !< convergence of lateral advective tracer fluxes
111 
112  integer :: ntr
113  type(tracer_type) :: temp
114  character(len=256) :: mesg ! Message for error messages.
115 
116  if (.not. associated(reg)) call tracer_registry_init(param_file, reg)
117 
118  if (reg%ntr>=max_fields_) then
119  write(mesg,'("Increase MAX_FIELDS_ in MOM_memory.h to at least ",I3," to allow for &
120  &all the tracers being registered via register_tracer.")') reg%ntr+1
121  call mom_error(fatal,"MOM register_tracer: "//mesg)
122  endif
123  reg%ntr = reg%ntr + 1
124  ntr = reg%ntr
125 
126  if (present(tr_desc_ptr)) then
127  reg%Tr(ntr)%vd => tr_desc_ptr
128  else
129  allocate(reg%Tr(ntr)%vd) ; reg%Tr(ntr)%vd = tr_desc
130  endif
131 
132  call query_vardesc(reg%Tr(ntr)%vd, name=reg%Tr(ntr)%name)
133 
134  if (reg%locked) call mom_error(fatal, &
135  "MOM register_tracer was called for variable "//trim(reg%Tr(ntr)%name)//&
136  " with a locked tracer registry.")
137 
138  reg%Tr(ntr)%t => tr1
139 
140  if (present(ad_x)) then ; if (associated(ad_x)) reg%Tr(ntr)%ad_x => ad_x ; endif
141  if (present(ad_y)) then ; if (associated(ad_y)) reg%Tr(ntr)%ad_y => ad_y ; endif
142  if (present(df_x)) then ; if (associated(df_x)) reg%Tr(ntr)%df_x => df_x ; endif
143  if (present(df_y)) then ; if (associated(df_y)) reg%Tr(ntr)%df_y => df_y ; endif
144  if (present(obc_inflow)) reg%Tr(ntr)%OBC_inflow_conc = obc_inflow
145  if (present(obc_in_u)) then ; if (associated(obc_in_u)) &
146  reg%Tr(ntr)%OBC_in_u => obc_in_u ; endif
147  if (present(obc_in_v)) then ; if (associated(obc_in_v)) &
148  reg%Tr(ntr)%OBC_in_v => obc_in_v ; endif
149  if (present(ad_2d_x)) then ; if (associated(ad_2d_x)) reg%Tr(ntr)%ad2d_x => ad_2d_x ; endif
150  if (present(ad_2d_y)) then ; if (associated(ad_2d_y)) reg%Tr(ntr)%ad2d_y => ad_2d_y ; endif
151  if (present(df_2d_x)) then ; if (associated(df_2d_x)) reg%Tr(ntr)%df2d_x => df_2d_x ; endif
152 
153  if (present(advection_xy)) then ; if (associated(advection_xy)) reg%Tr(ntr)%advection_xy => advection_xy ; endif
154 
155 end subroutine register_tracer
156 
157 
158 !> This subroutine locks the tracer registry to prevent the addition of more
159 !! tracers. After locked=.true., can then register common diagnostics.
160 subroutine lock_tracer_registry(Reg)
161  type(tracer_registry_type), pointer :: Reg !< pointer to the tracer registry
162 
163  if (.not. associated(reg)) call mom_error(warning, &
164  "lock_tracer_registry called with an unassocaited registry.")
165 
166  reg%locked = .true.
167 
168 end subroutine lock_tracer_registry
169 
170 
171 !> This subroutine adds open boundary condition concentrations for a tracer that
172 !! has previously been registered by a call to register_tracer.
173 subroutine add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
174  character(len=*), intent(in) :: name !< tracer name for which the diagnostic points
175  type(tracer_registry_type), pointer :: Reg !< pointer to the tracer registry
176  real, intent(in), optional :: OBC_inflow !< tracer value for all inflows via the OBC
177  !! for which OBC_in_u or OBC_in_v are
178  !! not specified (same units as tracer CONC)
179  real, pointer, dimension(:,:,:), optional :: OBC_in_u !< tracer at inflows through u-face of tracer cells
180  !! (same units as tracer CONC)
181  real, pointer, dimension(:,:,:), optional :: OBC_in_v !< tracer at inflows through v-face of tracer cells
182  !! (same units as tracer CONC)
183 
184  integer :: m
185 
186  if (.not. associated(reg)) call mom_error(fatal, "add_tracer_OBC_values :"// &
187  "register_tracer must be called before add_tracer_OBC_values")
188 
189  do m=1,reg%ntr ; if (reg%Tr(m)%name == trim(name)) exit ; enddo
190 
191  if (m <= reg%ntr) then
192  if (present(obc_inflow)) reg%Tr(m)%OBC_inflow_conc = obc_inflow
193  if (present(obc_in_u)) then ; if (associated(obc_in_u)) &
194  reg%Tr(m)%OBC_in_u => obc_in_u ; endif
195  if (present(obc_in_v)) then ; if (associated(obc_in_v)) &
196  reg%Tr(m)%OBC_in_v => obc_in_v ; endif
197  else
198  call mom_error(fatal, "MOM_tracer: register_tracer must be called for "//&
199  trim(name)//" before add_tracer_OBC_values is called for it.")
200  endif
201 
202 end subroutine add_tracer_obc_values
203 
204 
205 !> This subroutine adds diagnostic arrays for a tracer that has
206 !! previously been registered by a call to register_tracer.
207 subroutine add_tracer_diagnostics(name, Reg, ad_x, ad_y, df_x, df_y, &
208  ad_2d_x, ad_2d_y, df_2d_x, df_2d_y,&
209  advection_xy)
210  character(len=*), intent(in) :: name !< name of the tracer for which the diagnostic points
211  type(tracer_registry_type), pointer :: Reg !< pointer to the tracer registry
212  real, dimension(:,:,:), pointer, optional :: ad_x !< diagnostic x-advective flux (CONC m3/s or CONC*kg/s)
213  real, dimension(:,:,:), pointer, optional :: ad_y !< diagnostic y-advective flux (CONC m3/s or CONC*kg/s)
214  real, dimension(:,:,:), pointer, optional :: df_x !< diagnostic x-diffusive flux (CONC m3/s or CONC*kg/s)
215  real, dimension(:,:,:), pointer, optional :: df_y !< diagnostic y-diffusive flux (CONC m3/s or CONC*kg/s)
216  real, dimension(:,:), pointer, optional :: ad_2d_x !< vert sum of diagnostic x-advect flux (CONC m3/s or CONC*kg/s)
217  real, dimension(:,:), pointer, optional :: ad_2d_y !< vert sum of diagnostic y-advect flux (CONC m3/s or CONC*kg/s)
218  real, dimension(:,:), pointer, optional :: df_2d_x !< vert sum of diagnostic x-diffuse flux (CONC m3/s or CONC*kg/s)
219  real, dimension(:,:), pointer, optional :: df_2d_y !< vert sum of diagnostic y-diffuse flux (CONC m3/s or CONC*kg/s)
220 
221  real, dimension(:,:,:), pointer, optional :: advection_xy !< convergence of lateral advective tracer fluxes
222 
223  integer :: m
224 
225  if (.not. associated(reg)) call mom_error(fatal, "add_tracer_diagnostics: "// &
226  "register_tracer must be called before add_tracer_diagnostics")
227 
228  do m=1,reg%ntr ; if (reg%Tr(m)%name == trim(name)) exit ; enddo
229 
230  if (m <= reg%ntr) then
231  if (present(ad_x)) then ; if (associated(ad_x)) reg%Tr(m)%ad_x => ad_x ; endif
232  if (present(ad_y)) then ; if (associated(ad_y)) reg%Tr(m)%ad_y => ad_y ; endif
233  if (present(df_x)) then ; if (associated(df_x)) reg%Tr(m)%df_x => df_x ; endif
234  if (present(df_y)) then ; if (associated(df_y)) reg%Tr(m)%df_y => df_y ; endif
235 
236  if (present(ad_2d_x)) then ; if (associated(ad_2d_x)) reg%Tr(m)%ad2d_x => ad_2d_x ; endif
237  if (present(ad_2d_y)) then ; if (associated(ad_2d_y)) reg%Tr(m)%ad2d_y => ad_2d_y ; endif
238  if (present(df_2d_x)) then ; if (associated(df_2d_x)) reg%Tr(m)%df2d_x => df_2d_x ; endif
239  if (present(df_2d_y)) then ; if (associated(df_2d_y)) reg%Tr(m)%df2d_y => df_2d_y ; endif
240 
241  if (present(advection_xy)) then ; if (associated(advection_xy)) reg%Tr(m)%advection_xy => advection_xy ; endif
242 
243  else
244 
245  call mom_error(fatal, "MOM_tracer: register_tracer must be called for "//&
246  trim(name)//" before add_tracer_diagnostics is called for it.")
247  endif
248 
249 end subroutine add_tracer_diagnostics
250 
251 !> This subroutine writes out chksums for tracers.
252 subroutine mom_tracer_chksum(mesg, Tr, ntr, G)
253  character(len=*), intent(in) :: mesg !< message that appears on the chksum lines
254  type(tracer_type), intent(in) :: Tr(:) !< array of all of registered tracers
255  integer, intent(in) :: ntr !< number of registered tracers
256  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
257 
258  integer :: is, ie, js, je, nz
259  integer :: i, j, k, m
260 
261  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
262  do m=1,ntr
263  call hchksum(tr(m)%t, mesg//trim(tr(m)%name), g%HI)
264  enddo
265 
266 end subroutine mom_tracer_chksum
267 
268 !> Calculates and prints the global inventory of all tracers in the registry.
269 subroutine mom_tracer_chkinv(mesg, G, h, Tr, ntr)
270  character(len=*), intent(in) :: mesg !< message that appears on the chksum lines
271  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
272  type(tracer_type), intent(in) :: Tr(:) !< array of all of registered tracers
273  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses
274  integer, intent(in) :: ntr !< number of registered tracers
275 
276  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: tr_inv !< Tracer inventory
277  real :: total_inv
278  integer :: is, ie, js, je, nz
279  integer :: i, j, k, m
280 
281  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
282  do m=1,ntr
283  do k=1,nz ; do j=js,je ; do i=is,ie
284  tr_inv(i,j,k) = tr(m)%t(i,j,k)*h(i,j,k)*g%areaT(i,j)*g%mask2dT(i,j)
285  enddo ; enddo ; enddo
286  total_inv = reproducing_sum(tr_inv, is, ie, js, je)
287  if (is_root_pe()) write(0,'(A,1X,A5,1X,ES25.16,1X,A)') "h-point: inventory", tr(m)%name, total_inv, mesg
288  enddo
289 
290 end subroutine mom_tracer_chkinv
291 
292 !> This routine include declares and sets the variable "version".
293 subroutine tracer_registry_init(param_file, Reg)
294  type(param_file_type), intent(in) :: param_file !< open file to parse for model parameters
295  type(tracer_registry_type), pointer :: Reg !< pointer to tracer registry
296 
297  integer, save :: init_calls = 0
298 
299 #include "version_variable.h"
300  character(len=40) :: mdl = "MOM_tracer_registry" ! This module's name.
301  character(len=256) :: mesg ! Message for error messages.
302 
303  if (.not.associated(reg)) then ; allocate(reg)
304  else ; return ; endif
305 
306  ! Read all relevant parameters and write them to the model log.
307  call log_version(param_file, mdl, version, "")
308 
309  init_calls = init_calls + 1
310  if (init_calls > 1) then
311  write(mesg,'("tracer_registry_init called ",I3, &
312  &" times with different registry pointers.")') init_calls
313  if (is_root_pe()) call mom_error(warning,"MOM_tracer"//mesg)
314  endif
315 
316 end subroutine tracer_registry_init
317 
318 
319 !> This routine closes the tracer registry module.
320 subroutine tracer_registry_end(Reg)
321  type(tracer_registry_type), pointer :: Reg
322  if (associated(reg)) deallocate(reg)
323 end subroutine tracer_registry_end
324 
325 end module mom_tracer_registry
subroutine, public tracer_registry_init(param_file, Reg)
This routine include declares and sets the variable "version".
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Provides the ocean grid type.
Definition: MOM_grid.F90:2
subroutine, public mom_tracer_chksum(mesg, Tr, ntr, G)
This subroutine writes out chksums for tracers.
This module contains I/O framework code.
Definition: MOM_io.F90:2
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Container for horizontal index ranges for data, computational and global domains. ...
subroutine, public register_tracer(tr1, tr_desc, param_file, HI, GV, Reg, tr_desc_ptr, ad_x, ad_y, df_x, df_y, OBC_inflow, OBC_in_u, OBC_in_v, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine registers a tracer to be advected and laterally diffused.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Type to carry basic tracer information.
logical function, public is_root_pe()
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public add_tracer_diagnostics(name, Reg, ad_x, ad_y, df_x, df_y, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine adds diagnostic arrays for a tracer that has previously been registered by a call to ...
subroutine, public mom_mesg(message, verb, all_print)
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
subroutine, public tracer_registry_end(Reg)
This routine closes the tracer registry module.
subroutine, public lock_tracer_registry(Reg)
This subroutine locks the tracer registry to prevent the addition of more tracers. After locked=.true., can then register common diagnostics.
subroutine, public mom_tracer_chkinv(mesg, G, h, Tr, ntr)
Calculates and prints the global inventory of all tracers in the registry.
subroutine, public query_vardesc(vd, name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
This routine queries vardesc.
Definition: MOM_io.F90:664
subroutine, public mom_error(level, message, all_print)