MOM6
mom_generic_tracer Module Reference

Detailed Description

Drives the generic version of tracers TOPAZ and CFC and other GFDL BGC components.

Data Types

type  mom_generic_tracer_cs
 Control structure for generic tracers. More...
 

Functions/Subroutines

logical function, public register_mom_generic_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Initializes the generic tracer packages and adds their tracers to the list Adds the tracers in the list of generic tracers to the set of MOM tracers (i.e., MOM-register them) Register these tracers for restart. More...
 
subroutine, public initialize_mom_generic_tracer (restart, day, G, GV, US, h, param_file, diag, OBC, CS, sponge_CSp, ALE_sponge_CSp)
 Initialize phase II: Initialize required variables for generic tracers There are some steps of initialization that cannot be done in register_MOM_generic_tracer This is the place and time to do them: Set the grid mask and initial time for all generic tracers. Diag_register them. Z_diag_register them. More...
 
subroutine, public mom_generic_tracer_column_physics (h_old, h_new, ea, eb, fluxes, Hml, dt, G, GV, CS, tv, optics, evap_CFL_limit, minimum_forcing_depth)
 Column physics for generic tracers. Get the coupler values for generic tracers that exchange with atmosphere Update generic tracer concentration fields from sources and sinks. Vertically diffuse generic tracer concentration fields. Update generic tracers from bottom and their bottom reservoir. More...
 
integer function, public mom_generic_tracer_stock (h, stocks, G, GV, CS, names, units, stock_index)
 This subroutine calculates mass-weighted integral on the PE either of all available tracer concentrations, or of a tracer that is being requested specifically, returning the number of stocks it has calculated. If the stock_index is present, only the stock corresponding to that coded index is returned. More...
 
integer function, public mom_generic_tracer_min_max (ind_start, got_minmax, gmin, gmax, xgmin, ygmin, zgmin, xgmax, ygmax, zgmax, G, CS, names, units)
 This subroutine find the global min and max of either of all available tracer concentrations, or of a tracer that is being requested specifically, returning the number of tracers it has gone through. More...
 
subroutine, public mom_generic_tracer_surface_state (sfc_state, h, G, CS)
 This subroutine calculates the surface state and sets coupler values for those generic tracers that have flux exchange with atmosphere. More...
 
subroutine, public mom_generic_flux_init (verbosity)
 
subroutine, public mom_generic_tracer_fluxes_accumulate (flux_tmp, weight)
 
subroutine, public mom_generic_tracer_get (name, member, array, CS)
 Copy the requested tracer into an array. More...
 
subroutine, public end_mom_generic_tracer (CS)
 This subroutine deallocates the memory owned by this module. More...
 

Variables

logical g_registered = .false.
 An state hidden in module data that is very much not allowed in MOM6.
 

Function/Subroutine Documentation

◆ end_mom_generic_tracer()

subroutine, public mom_generic_tracer::end_mom_generic_tracer ( type(mom_generic_tracer_cs), pointer  CS)

This subroutine deallocates the memory owned by this module.

Parameters
csPointer to the control structure for this module.

Definition at line 822 of file MOM_generic_tracer.F90.

822  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
823 
824  call generic_tracer_end()
825 
826  if (associated(cs)) then
827  deallocate(cs)
828  endif

◆ initialize_mom_generic_tracer()

subroutine, public mom_generic_tracer::initialize_mom_generic_tracer ( logical, intent(in)  restart,
type(time_type), intent(in), target  day,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(mom_generic_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp,
type(ale_sponge_cs), pointer  ALE_sponge_CSp 
)

Initialize phase II: Initialize required variables for generic tracers There are some steps of initialization that cannot be done in register_MOM_generic_tracer This is the place and time to do them: Set the grid mask and initial time for all generic tracers. Diag_register them. Z_diag_register them.

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

Parameters
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in,out]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]param_fileA structure to parse for run-time parameters
[in]diagRegulates diagnostic output.
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used.
csPointer to the control structure for this module.
sponge_cspPointer to the control structure for the sponges.
ale_sponge_cspPointer to the control structure for the ALE sponges.

Definition at line 232 of file MOM_generic_tracer.F90.

232  logical, intent(in) :: restart !< .true. if the fields have already been
233  !! read from a restart file.
234  type(time_type), target, intent(in) :: day !< Time of the start of the run.
235  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
236  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
237  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
238  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
239  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
240  type(diag_ctrl), target, intent(in) :: diag !< Regulates diagnostic output.
241  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies whether,
242  !! where, and what open boundary conditions are used.
243  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
244  type(sponge_cs), pointer :: sponge_csp !< Pointer to the control structure for the sponges.
245  type(ale_sponge_cs), pointer :: ale_sponge_csp !< Pointer to the control structure for the
246  !! ALE sponges.
247 
248  character(len=128), parameter :: sub_name = 'initialize_MOM_generic_tracer'
249  logical :: ok
250  integer :: i, j, k, isc, iec, jsc, jec, nk
251  type(g_tracer_type), pointer :: g_tracer,g_tracer_next
252  character(len=fm_string_len) :: g_tracer_name
253  real, dimension(:,:,:,:), pointer :: tr_field
254  real, dimension(:,:,:), pointer :: tr_ptr
255  real, dimension(G%isd:G%ied, G%jsd:G%jed,1:G%ke) :: grid_tmask
256  integer, dimension(G%isd:G%ied, G%jsd:G%jed) :: grid_kmt
257 
258  !! 2010/02/04 Add code to re-initialize Generic Tracers if needed during a model simulation
259  !! By default, restart cpio should not contain a Generic Tracer IC file and step below will be skipped.
260  !! Ideally, the generic tracer IC file should have the tracers on Z levels.
261 
262  isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec ; nk = g%ke
263 
264  cs%diag=>diag
265  !Get the tracer list
266  if (.NOT. associated(cs%g_tracer_list)) call mom_error(fatal, trim(sub_name)//&
267  ": No tracer in the list.")
268  !For each tracer name get its fields
269  g_tracer=>cs%g_tracer_list
270 
271  do
272  if (index(cs%IC_file, '_NULL_') /= 0) then
273  call mom_error(warning, "The name of the IC_file "//trim(cs%IC_file)//&
274  " indicates no MOM initialization was asked for the generic tracers."//&
275  "Bypassing the MOM initialization of ALL generic tracers!")
276  exit
277  endif
278  call g_tracer_get_alias(g_tracer,g_tracer_name)
279  call g_tracer_get_pointer(g_tracer,g_tracer_name,'field',tr_field)
280  tr_ptr => tr_field(:,:,:,1)
281 
282  if (.not.restart .or. (cs%tracers_may_reinit .and. &
283  .not.query_initialized(tr_ptr, g_tracer_name, cs%restart_CSp))) then
284 
285  if (g_tracer%requires_src_info ) then
286  call mom_error(note,"initialize_MOM_generic_tracer: "//&
287  "initializing generic tracer "//trim(g_tracer_name)//&
288  " using MOM_initialize_tracer_from_Z ")
289 
290  call mom_initialize_tracer_from_z(h, tr_ptr, g, gv, us, param_file, &
291  src_file = g_tracer%src_file, &
292  src_var_nam = g_tracer%src_var_name, &
293  src_var_unit_conversion = g_tracer%src_var_unit_conversion,&
294  src_var_record = g_tracer%src_var_record, &
295  src_var_gridspec = g_tracer%src_var_gridspec )
296 
297  !Check/apply the bounds for each g_tracer
298  do k=1,nk ; do j=jsc,jec ; do i=isc,iec
299  if (tr_ptr(i,j,k) /= cs%tracer_land_val) then
300  if (tr_ptr(i,j,k) < g_tracer%src_var_valid_min) tr_ptr(i,j,k) = g_tracer%src_var_valid_min
301  !Jasmin does not want to apply the maximum for now
302  !if (tr_ptr(i,j,k) > g_tracer%src_var_valid_max) tr_ptr(i,j,k) = g_tracer%src_var_valid_max
303  endif
304  enddo ; enddo ; enddo
305 
306  !jgj: Reset CASED to 0 below K=1
307  if ( (trim(g_tracer_name) == 'cased') .or. (trim(g_tracer_name) == 'ca13csed') ) then
308  do k=2,nk ; do j=jsc,jec ; do i=isc,iec
309  if (tr_ptr(i,j,k) /= cs%tracer_land_val) then
310  tr_ptr(i,j,k) = 0.0
311  endif
312  enddo ; enddo ; enddo
313  endif
314  elseif(.not. g_tracer%requires_restart) then
315  !Do nothing for this tracer, it is initialized by the tracer package
316  call mom_error(note,"initialize_MOM_generic_tracer: "//&
317  "skip initialization of generic tracer "//trim(g_tracer_name))
318  else !Do it old way if the tracer is not registered to start from a specific source file.
319  !This path should be deprecated if all generic tracers are required to start from specified sources.
320  if (len_trim(cs%IC_file) > 0) then
321  ! Read the tracer concentrations from a netcdf file.
322  if (.not.file_exists(cs%IC_file)) call mom_error(fatal, &
323  "initialize_MOM_Generic_tracer: Unable to open "//cs%IC_file)
324  if (cs%Z_IC_file) then
325  ok = tracer_z_init(tr_ptr, h, cs%IC_file, g_tracer_name, g, us)
326  if (.not.ok) then
327  ok = tracer_z_init(tr_ptr, h, cs%IC_file, trim(g_tracer_name), g, us)
328  if (.not.ok) call mom_error(fatal,"initialize_MOM_Generic_tracer: "//&
329  "Unable to read "//trim(g_tracer_name)//" from "//&
330  trim(cs%IC_file)//".")
331  endif
332  call mom_error(note,"initialize_MOM_generic_tracer: "//&
333  "initialized generic tracer "//trim(g_tracer_name)//&
334  " using Generic Tracer File on Z: "//cs%IC_file)
335  else
336  ! native grid
337  call mom_error(note,"initialize_MOM_generic_tracer: "//&
338  "Using Generic Tracer IC file on native grid "//trim(cs%IC_file)//&
339  " for tracer "//trim(g_tracer_name))
340  call mom_read_data(cs%IC_file, trim(g_tracer_name), tr_ptr, g%Domain)
341  endif
342  else
343  call mom_error(fatal,"initialize_MOM_generic_tracer: "//&
344  "check Generic Tracer IC filename "//trim(cs%IC_file)//&
345  " for tracer "//trim(g_tracer_name))
346  endif
347 
348  endif
349  endif
350 
351  !traverse the linked list till hit NULL
352  call g_tracer_get_next(g_tracer, g_tracer_next)
353  if (.NOT. associated(g_tracer_next)) exit
354  g_tracer=>g_tracer_next
355  enddo
356  !! end section to re-initialize generic tracers
357 
358 
359  !Now we can reset the grid mask, axes and time to their true values
360  !Note that grid_tmask must be set correctly on the data domain boundary
361  !so that coast mask can be deduced from it.
362  grid_tmask(:,:,:) = 0.0
363  grid_kmt(:,:) = 0
364  do j = g%jsd, g%jed ; do i = g%isd, g%ied
365  if (g%mask2dT(i,j) > 0) then
366  grid_tmask(i,j,:) = 1.0
367  grid_kmt(i,j) = g%ke ! Tell the code that a layer thicker than 1m is the bottom layer.
368  endif
369  enddo ; enddo
370  call g_tracer_set_common(g%isc,g%iec,g%jsc,g%jec,g%isd,g%ied,g%jsd,g%jed,&
371  gv%ke,1,cs%diag%axesTL%handles,grid_tmask,grid_kmt,day)
372 
373  ! Register generic tracer modules diagnostics
374 
375 #ifdef _USE_MOM6_DIAG
376  call g_tracer_set_csdiag(cs%diag)
377 #endif
378  call generic_tracer_register_diag()
379 #ifdef _USE_MOM6_DIAG
380  call g_tracer_set_csdiag(cs%diag)
381 #endif
382 
383  cs%H_to_m = gv%H_to_m
384 

◆ mom_generic_flux_init()

subroutine, public mom_generic_tracer::mom_generic_flux_init ( integer, intent(in), optional  verbosity)
Parameters
[in]verbosityA 0-9 integer indicating a level of verbosity.

Definition at line 763 of file MOM_generic_tracer.F90.

763  integer, optional, intent(in) :: verbosity !< A 0-9 integer indicating a level of verbosity.
764 
765  integer :: ind
766  character(len=fm_string_len) :: g_tracer_name,longname, package,units,old_package,file_in,file_out
767  real :: const_init_value
768  character(len=128), parameter :: sub_name = 'MOM_generic_flux_init'
769  type(g_tracer_type), pointer :: g_tracer_list,g_tracer,g_tracer_next
770 
771  if (.not. g_registered) then
772  call generic_tracer_register()
773  g_registered = .true.
774  endif
775 
776  call generic_tracer_get_list(g_tracer_list)
777  if (.NOT. associated(g_tracer_list)) then
778  call mom_error(warning, trim(sub_name)// ": No generic tracer in the list.")
779  return
780  endif
781 
782  g_tracer=>g_tracer_list
783  do
784 
785  call g_tracer_flux_init(g_tracer) !, verbosity=verbosity) !### Add this after ocean shared is updated.
786 
787  ! traverse the linked list till hit NULL
788  call g_tracer_get_next(g_tracer, g_tracer_next)
789  if (.NOT. associated(g_tracer_next)) exit
790  g_tracer=>g_tracer_next
791 
792  enddo
793 

◆ mom_generic_tracer_column_physics()

subroutine, public mom_generic_tracer::mom_generic_tracer_column_physics ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_old,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_new,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  ea,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  eb,
type(forcing), intent(in)  fluxes,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(in)  Hml,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(mom_generic_tracer_cs), pointer  CS,
type(thermo_var_ptrs), intent(in)  tv,
type(optics_type), intent(in)  optics,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

Column physics for generic tracers. Get the coupler values for generic tracers that exchange with atmosphere Update generic tracer concentration fields from sources and sinks. Vertically diffuse generic tracer concentration fields. Update generic tracers from bottom and their bottom reservoir.

This subroutine applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this file. CFCs are relatively simple, as they are passive tracers. with only a surface flux as a source.

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]eaThe amount of fluid entrained from the layer
[in]ebThe amount of fluid entrained from the layer
[in]fluxesA structure containing pointers to thermodynamic and tracer forcing fields.
[in]hmlMixed layer depth [Z ~> m]
[in]dtThe amount of time covered by this call [s]
csPointer to the control structure for this module.
[in]tvA structure pointing to various thermodynamic variables
[in]opticsThe structure containing optical properties.
[in]evap_cfl_limitLimits how much water can be fluxed out of the top layer Stored previously in diabatic CS.
[in]minimum_forcing_depthThe smallest depth over which fluxes can be applied [H ~> m or kg m-2]

Definition at line 399 of file MOM_generic_tracer.F90.

399  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
400  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
401  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
402  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
403  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
404  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
405  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
406  intent(in) :: ea !< The amount of fluid entrained from the layer
407  !! above during this call [H ~> m or kg m-2].
408  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
409  intent(in) :: eb !< The amount of fluid entrained from the layer
410  !! below during this call [H ~> m or kg m-2].
411  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
412  !! and tracer forcing fields.
413  real, dimension(SZI_(G),SZJ_(G)), intent(in) :: hml !< Mixed layer depth [Z ~> m]
414  real, intent(in) :: dt !< The amount of time covered by this call [s]
415  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
416  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
417  type(optics_type), intent(in) :: optics !< The structure containing optical properties.
418  real, optional, intent(in) :: evap_cfl_limit !< Limits how much water can be fluxed out of
419  !! the top layer Stored previously in diabatic CS.
420  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which fluxes
421  !! can be applied [H ~> m or kg m-2]
422  ! Stored previously in diabatic CS.
423  ! The arguments to this subroutine are redundant in that
424  ! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
425 
426  ! Local variables
427  character(len=128), parameter :: sub_name = 'MOM_generic_tracer_column_physics'
428 
429  type(g_tracer_type), pointer :: g_tracer, g_tracer_next
430  character(len=fm_string_len) :: g_tracer_name
431  real, dimension(:,:), pointer :: stf_array,trunoff_array,runoff_tracer_flux_array
432 
433  real :: surface_field(szi_(g),szj_(g))
434  real :: dz_ml(szi_(g),szj_(g)) ! The mixed layer depth in the MKS units used for generic tracers [m]
435  real :: sosga
436 
437  real, dimension(G%isd:G%ied,G%jsd:G%jed,G%ke) :: rho_dzt, dzt
438  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work
439  integer :: i, j, k, isc, iec, jsc, jec, nk
440 
441  isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec ; nk = g%ke
442 
443  !Get the tracer list
444  if (.NOT. associated(cs%g_tracer_list)) call mom_error(fatal,&
445  trim(sub_name)//": No tracer in the list.")
446 
447 #ifdef _USE_MOM6_DIAG
448  call g_tracer_set_csdiag(cs%diag)
449 #endif
450 
451  !
452  !Extract the tracer surface fields from coupler and update tracer fields from sources
453  !
454  !call generic_tracer_coupler_get(fluxes%tr_fluxes)
455  !Niki: This is moved out to ocean_model_MOM.F90 because if dt_therm>dt_cpld we need to average
456  ! the fluxes without coming into this subroutine.
457  ! MOM5 has to modified to conform.
458 
459  !
460  !Add contribution of river to surface flux
461  !
462  g_tracer=>cs%g_tracer_list
463  do
464  if (_allocated(g_tracer%trunoff)) then
465  call g_tracer_get_alias(g_tracer,g_tracer_name)
466  call g_tracer_get_pointer(g_tracer,g_tracer_name,'stf', stf_array)
467  call g_tracer_get_pointer(g_tracer,g_tracer_name,'trunoff',trunoff_array)
468  call g_tracer_get_pointer(g_tracer,g_tracer_name,'runoff_tracer_flux',runoff_tracer_flux_array)
469  !nnz: Why is fluxes%river = 0?
470  runoff_tracer_flux_array(:,:) = trunoff_array(:,:) * &
471  g%US%RZ_T_to_kg_m2s*fluxes%lrunoff(:,:)
472  stf_array = stf_array + runoff_tracer_flux_array
473  endif
474 
475  !traverse the linked list till hit NULL
476  call g_tracer_get_next(g_tracer, g_tracer_next)
477  if (.NOT. associated(g_tracer_next)) exit
478  g_tracer=>g_tracer_next
479 
480  enddo
481 
482  !
483  !Prepare input arrays for source update
484  !
485 
486  rho_dzt(:,:,:) = gv%H_to_kg_m2 * gv%Angstrom_H
487  do k = 1, nk ; do j = jsc, jec ; do i = isc, iec !{
488  rho_dzt(i,j,k) = gv%H_to_kg_m2 * h_old(i,j,k)
489  enddo ; enddo ; enddo !}
490 
491  dzt(:,:,:) = 1.0
492  do k = 1, nk ; do j = jsc, jec ; do i = isc, iec !{
493  dzt(i,j,k) = gv%H_to_m * h_old(i,j,k)
494  enddo ; enddo ; enddo !}
495  dz_ml(:,:) = 0.0
496  do j=jsc,jec ; do i=isc,iec
497  surface_field(i,j) = tv%S(i,j,1)
498  dz_ml(i,j) = g%US%Z_to_m * hml(i,j)
499  enddo ; enddo
500  sosga = global_area_mean(surface_field, g)
501 
502  !
503  !Calculate tendencies (i.e., field changes at dt) from the sources / sinks
504  !
505  if ((g%US%L_to_m == 1.0) .and. (g%US%RZ_to_kg_m2 == 1.0) .and. (g%US%s_to_T == 1.0)) then
506  ! Avoid unnecessary copies when no unit conversion is needed.
507  call generic_tracer_source(tv%T, tv%S, rho_dzt, dzt, dz_ml, g%isd, g%jsd, 1, dt, &
508  g%areaT, get_diag_time_end(cs%diag), &
509  optics%nbands, optics%max_wavelength_band, optics%sw_pen_band, optics%opacity_band, &
510  internal_heat=tv%internal_heat, frunoff=fluxes%frunoff, sosga=sosga)
511  else
512  call generic_tracer_source(tv%T, tv%S, rho_dzt, dzt, dz_ml, g%isd, g%jsd, 1, dt, &
513  g%US%L_to_m**2*g%areaT(:,:), get_diag_time_end(cs%diag), &
514  optics%nbands, optics%max_wavelength_band, optics%sw_pen_band, optics%opacity_band, &
515  internal_heat=g%US%RZ_to_kg_m2*tv%internal_heat(:,:), &
516  frunoff=g%US%RZ_T_to_kg_m2s*fluxes%frunoff(:,:), sosga=sosga)
517  endif
518 
519  ! This uses applyTracerBoundaryFluxesInOut to handle the change in tracer due to freshwater fluxes
520  ! usually in ALE mode
521  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
522  g_tracer=>cs%g_tracer_list
523  do
524  if (g_tracer_is_prog(g_tracer)) then
525  do k=1,nk ;do j=jsc,jec ; do i=isc,iec
526  h_work(i,j,k) = h_old(i,j,k)
527  enddo ; enddo ; enddo
528  call applytracerboundaryfluxesinout(g, gv, g_tracer%field(:,:,:,1), g%US%s_to_T*dt, &
529  fluxes, h_work, evap_cfl_limit, minimum_forcing_depth)
530  endif
531 
532  !traverse the linked list till hit NULL
533  call g_tracer_get_next(g_tracer, g_tracer_next)
534  if (.NOT. associated(g_tracer_next)) exit
535  g_tracer=>g_tracer_next
536  enddo
537  endif
538 
539  !
540  !Update Tr(n)%field from explicit vertical diffusion
541  !
542  ! Use a tridiagonal solver to determine the concentrations after the
543  ! surface source is applied and diapycnal advection and diffusion occurs.
544  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
545  ! Last arg is tau which is always 1 for MOM6
546  call generic_tracer_vertdiff_g(h_work, ea, eb, dt, gv%kg_m2_to_H, gv%m_to_H, 1)
547  else
548  ! Last arg is tau which is always 1 for MOM6
549  call generic_tracer_vertdiff_g(h_old, ea, eb, dt, gv%kg_m2_to_H, gv%m_to_H, 1)
550  endif
551 
552  ! Update bottom fields after vertical processes
553 
554  ! Second arg is tau which is always 1 for MOM6
555  call generic_tracer_update_from_bottom(dt, 1, get_diag_time_end(cs%diag))
556 
557  !Output diagnostics via diag_manager for all generic tracers and their fluxes
558  call g_tracer_send_diag(cs%g_tracer_list, get_diag_time_end(cs%diag), tau=1)
559 #ifdef _USE_MOM6_DIAG
560  call g_tracer_set_csdiag(cs%diag)
561 #endif
562 

◆ mom_generic_tracer_fluxes_accumulate()

subroutine, public mom_generic_tracer::mom_generic_tracer_fluxes_accumulate ( type(forcing), intent(in)  flux_tmp,
real, intent(in)  weight 
)
Parameters
[in]flux_tmpA structure containing pointers to thermodynamic and tracer forcing fields.
[in]weightA weight for accumulating this flux

Definition at line 797 of file MOM_generic_tracer.F90.

797  type(forcing), intent(in) :: flux_tmp !< A structure containing pointers to
798  !! thermodynamic and tracer forcing fields.
799  real, intent(in) :: weight !< A weight for accumulating this flux
800 
801  call generic_tracer_coupler_accumulate(flux_tmp%tr_fluxes, weight)
802 

◆ mom_generic_tracer_get()

subroutine, public mom_generic_tracer::mom_generic_tracer_get ( character(len=*), intent(in)  name,
character(len=*), intent(in)  member,
real, dimension(:,:,:), intent(out)  array,
type(mom_generic_tracer_cs), pointer  CS 
)

Copy the requested tracer into an array.

Parameters
[in]nameName of requested tracer.
[in]memberThe tracer element to return.
[out]arrayArray filled by this routine.
csPointer to the control structure for this module.

Definition at line 807 of file MOM_generic_tracer.F90.

807  character(len=*), intent(in) :: name !< Name of requested tracer.
808  character(len=*), intent(in) :: member !< The tracer element to return.
809  real, dimension(:,:,:), intent(out) :: array !< Array filled by this routine.
810  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
811 
812  real, dimension(:,:,:), pointer :: array_ptr
813  character(len=128), parameter :: sub_name = 'MOM_generic_tracer_get'
814 
815  call g_tracer_get_pointer(cs%g_tracer_list,name,member,array_ptr)
816  array(:,:,:) = array_ptr(:,:,:)
817 

◆ mom_generic_tracer_min_max()

integer function, public mom_generic_tracer::mom_generic_tracer_min_max ( integer, intent(in)  ind_start,
logical, dimension(:), intent(out)  got_minmax,
real, dimension(:), intent(out)  gmin,
real, dimension(:), intent(out)  gmax,
real, dimension(:), intent(out)  xgmin,
real, dimension(:), intent(out)  ygmin,
real, dimension(:), intent(out)  zgmin,
real, dimension(:), intent(out)  xgmax,
real, dimension(:), intent(out)  ygmax,
real, dimension(:), intent(out)  zgmax,
type(ocean_grid_type), intent(in)  G,
type(mom_generic_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units 
)

This subroutine find the global min and max of either of all available tracer concentrations, or of a tracer that is being requested specifically, returning the number of tracers it has gone through.

Parameters
[in]ind_startThe index of the tracer to start with
[out]got_minmaxIndicates whether the global min and max are found for each tracer
[out]gminGlobal minimum of each tracer, in kg times concentration units.
[out]gmaxGlobal maximum of each tracer, in kg times concentration units.
[out]xgminThe x-position of the global minimum
[out]ygminThe y-position of the global minimum
[out]zgminThe z-position of the global minimum
[out]xgmaxThe x-position of the global maximum
[out]ygmaxThe y-position of the global maximum
[out]zgmaxThe z-position of the global maximum
[in]gThe ocean's grid structure
csPointer to the control structure for this module.
[out]namesThe names of the stocks calculated.
[out]unitsThe units of the stocks calculated.
Returns
Return value, the number of tracers done here.

Definition at line 636 of file MOM_generic_tracer.F90.

636  use mpp_utilities_mod, only: mpp_array_global_min_max
637  integer, intent(in) :: ind_start !< The index of the tracer to start with
638  logical, dimension(:), intent(out) :: got_minmax !< Indicates whether the global min and
639  !! max are found for each tracer
640  real, dimension(:), intent(out) :: gmin !< Global minimum of each tracer, in kg
641  !! times concentration units.
642  real, dimension(:), intent(out) :: gmax !< Global maximum of each tracer, in kg
643  !! times concentration units.
644  real, dimension(:), intent(out) :: xgmin !< The x-position of the global minimum
645  real, dimension(:), intent(out) :: ygmin !< The y-position of the global minimum
646  real, dimension(:), intent(out) :: zgmin !< The z-position of the global minimum
647  real, dimension(:), intent(out) :: xgmax !< The x-position of the global maximum
648  real, dimension(:), intent(out) :: ygmax !< The y-position of the global maximum
649  real, dimension(:), intent(out) :: zgmax !< The z-position of the global maximum
650  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
651  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
652  character(len=*), dimension(:), intent(out) :: names !< The names of the stocks calculated.
653  character(len=*), dimension(:), intent(out) :: units !< The units of the stocks calculated.
654  integer :: mom_generic_tracer_min_max !< Return value, the
655  !! number of tracers done here.
656 
657 ! Local variables
658  type(g_tracer_type), pointer :: g_tracer, g_tracer_next
659  real, dimension(:,:,:,:), pointer :: tr_field
660  real, dimension(:,:,:), pointer :: tr_ptr
661  character(len=128), parameter :: sub_name = 'MOM_generic_tracer_min_max'
662 
663  real, dimension(:,:,:),pointer :: grid_tmask
664  integer :: isc,iec,jsc,jec,isd,ied,jsd,jed,nk,ntau
665 
666  integer :: i, j, k, is, ie, js, je, nz, m
667  real, allocatable, dimension(:) :: geo_z
668 
669  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
670 
671  mom_generic_tracer_min_max = 0
672  if (.not.associated(cs)) return
673 
674  if (.NOT. associated(cs%g_tracer_list)) return ! No stocks.
675 
676 
677  call g_tracer_get_common(isc,iec,jsc,jec,isd,ied,jsd,jed,nk,ntau,grid_tmask=grid_tmask)
678 
679  ! Because the use of a simple z-coordinate can not be assumed, simply
680  ! use the layer index as the vertical label.
681  allocate(geo_z(nk))
682  do k=1,nk ; geo_z(k) = real(k) ; enddo
683 
684  m=ind_start ; g_tracer=>cs%g_tracer_list
685  do
686  call g_tracer_get_alias(g_tracer,names(m))
687  call g_tracer_get_values(g_tracer,names(m),'units',units(m))
688  units(m) = trim(units(m))//" kg"
689  call g_tracer_get_pointer(g_tracer,names(m),'field',tr_field)
690 
691  gmin(m) = -1.0
692  gmax(m) = -1.0
693 
694  tr_ptr => tr_field(:,:,:,1)
695 
696  call mpp_array_global_min_max(tr_ptr, grid_tmask,isd,jsd,isc,iec,jsc,jec,nk , gmin(m), gmax(m), &
697  g%geoLonT,g%geoLatT,geo_z,xgmin(m), ygmin(m), zgmin(m), &
698  xgmax(m), ygmax(m), zgmax(m))
699 
700  got_minmax(m) = .true.
701 
702  !traverse the linked list till hit NULL
703  call g_tracer_get_next(g_tracer, g_tracer_next)
704  if (.NOT. associated(g_tracer_next)) exit
705  g_tracer=>g_tracer_next
706  m = m+1
707  enddo
708 
709  mom_generic_tracer_min_max = m
710 

◆ mom_generic_tracer_stock()

integer function, public mom_generic_tracer::mom_generic_tracer_stock ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
real, dimension(:), intent(out)  stocks,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(mom_generic_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units,
integer, intent(in), optional  stock_index 
)

This subroutine calculates mass-weighted integral on the PE either of all available tracer concentrations, or of a tracer that is being requested specifically, returning the number of stocks it has calculated. If the stock_index is present, only the stock corresponding to that coded index is returned.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses [H ~> m or kg m-2]
[out]stocksThe mass-weighted integrated amount of each tracer, in kg times concentration units [kg conc].
csPointer to the control structure for this module.
[out]namesThe names of the stocks calculated.
[out]unitsThe units of the stocks calculated.
[in]stock_indexThe coded index of a specific stock being sought.
Returns
Return value, the number of stocks calculated here.

Definition at line 571 of file MOM_generic_tracer.F90.

571  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
572  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
573  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
574  real, dimension(:), intent(out) :: stocks !< The mass-weighted integrated amount of each
575  !! tracer, in kg times concentration units [kg conc].
576  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
577  character(len=*), dimension(:), intent(out) :: names !< The names of the stocks calculated.
578  character(len=*), dimension(:), intent(out) :: units !< The units of the stocks calculated.
579  integer, optional, intent(in) :: stock_index !< The coded index of a specific stock
580  !! being sought.
581  integer :: mom_generic_tracer_stock !< Return value, the
582  !! number of stocks calculated here.
583 
584 ! Local variables
585  type(g_tracer_type), pointer :: g_tracer, g_tracer_next
586  real, dimension(:,:,:,:), pointer :: tr_field
587  real, dimension(:,:,:), pointer :: tr_ptr
588  character(len=128), parameter :: sub_name = 'MOM_generic_tracer_stock'
589 
590  integer :: i, j, k, is, ie, js, je, nz, m
591  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
592 
593  mom_generic_tracer_stock = 0
594  if (.not.associated(cs)) return
595 
596  if (present(stock_index)) then ; if (stock_index > 0) then
597  ! Check whether this stock is available from this routine.
598 
599  ! No stocks from this routine are being checked yet. Return 0.
600  return
601  endif ; endif
602 
603  if (.NOT. associated(cs%g_tracer_list)) return ! No stocks.
604 
605  m=1 ; g_tracer=>cs%g_tracer_list
606  do
607  call g_tracer_get_alias(g_tracer,names(m))
608  call g_tracer_get_values(g_tracer,names(m),'units',units(m))
609  units(m) = trim(units(m))//" kg"
610  call g_tracer_get_pointer(g_tracer,names(m),'field',tr_field)
611 
612  stocks(m) = 0.0
613  tr_ptr => tr_field(:,:,:,1)
614  do k=1,nz ; do j=js,je ; do i=is,ie
615  stocks(m) = stocks(m) + tr_ptr(i,j,k) * &
616  (g%mask2dT(i,j) * g%US%L_to_m**2*g%areaT(i,j) * h(i,j,k))
617  enddo ; enddo ; enddo
618  stocks(m) = gv%H_to_kg_m2 * stocks(m)
619 
620  !traverse the linked list till hit NULL
621  call g_tracer_get_next(g_tracer, g_tracer_next)
622  if (.NOT. associated(g_tracer_next)) exit
623  g_tracer=>g_tracer_next
624  m = m+1
625  enddo
626 
627  mom_generic_tracer_stock = m
628 

◆ mom_generic_tracer_surface_state()

subroutine, public mom_generic_tracer::mom_generic_tracer_surface_state ( type(surface), intent(inout)  sfc_state,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(mom_generic_tracer_cs), pointer  CS 
)

This subroutine calculates the surface state and sets coupler values for those generic tracers that have flux exchange with atmosphere.

This subroutine sets up the fields that the coupler needs to calculate the CFC fluxes between the ocean and atmosphere.

Parameters
[in]gThe ocean's grid structure
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in]hLayer thicknesses [H ~> m or kg m-2]
csPointer to the control structure for this module.

Definition at line 720 of file MOM_generic_tracer.F90.

720  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
721  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
722  !! describe the surface state of the ocean.
723  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
724  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module.
725 
726 ! Local variables
727  real :: sosga
728 
729  character(len=128), parameter :: sub_name = 'MOM_generic_tracer_surface_state'
730  real, dimension(G%isd:G%ied,G%jsd:G%jed,1:G%ke,1) :: rho0
731  real, dimension(G%isd:G%ied,G%jsd:G%jed,1:G%ke) :: dzt
732  type(g_tracer_type), pointer :: g_tracer
733 
734  !Set coupler values
735  !nnz: fake rho0
736  rho0=1.0
737 
738  dzt(:,:,:) = cs%H_to_m * h(:,:,:)
739 
740  sosga = global_area_mean(sfc_state%SSS, g)
741 
742  call generic_tracer_coupler_set(sfc_state%tr_fields,&
743  st=sfc_state%SST,&
744  ss=sfc_state%SSS,&
745  rho=rho0,& !nnz: required for MOM5 and previous versions.
746  ilb=g%isd, jlb=g%jsd,&
747  dzt=dzt,& !This is needed for the Mocsy method of carbonate system vars
748  tau=1,sosga=sosga,model_time=get_diag_time_end(cs%diag))
749 
750  !Output diagnostics via diag_manager for all tracers in this module
751 ! if (.NOT. associated(CS%g_tracer_list)) call MOM_error(FATAL, trim(sub_name)//&
752 ! "No tracer in the list.")
753 ! call g_tracer_send_diag(CS%g_tracer_list, get_diag_time_end(CS%diag), tau=1)
754  !Niki: The problem with calling diagnostic outputs here is that this subroutine is called every dt_cpld
755  ! hence if dt_therm > dt_cpld we get output (and contribution to the mean) at times that tracers
756  ! had not been updated.
757  ! Moving this to the end of column physics subrotuine fixes this issue.
758 

◆ register_mom_generic_tracer()

logical function, public mom_generic_tracer::register_mom_generic_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(mom_generic_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

Initializes the generic tracer packages and adds their tracers to the list Adds the tracers in the list of generic tracers to the set of MOM tracers (i.e., MOM-register them) Register these tracers for restart.

Parameters
[in]hiHorizontal index ranges
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters
csPointer to the control structure for this module
tr_regPointer to the control structure for the tracer advection and diffusion module.
restart_csPointer to the restart control structure.

Definition at line 98 of file MOM_generic_tracer.F90.

98  type(hor_index_type), intent(in) :: hi !< Horizontal index ranges
99  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
100  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
101  type(mom_generic_tracer_cs), pointer :: cs !< Pointer to the control structure for this module
102  type(tracer_registry_type), pointer :: tr_reg !< Pointer to the control structure for the tracer
103  !! advection and diffusion module.
104  type(mom_restart_cs), pointer :: restart_cs !< Pointer to the restart control structure.
105 
106 ! Local variables
107  logical :: register_mom_generic_tracer
108 
109  character(len=128), parameter :: sub_name = 'register_MOM_generic_tracer'
110  character(len=200) :: inputdir ! The directory where NetCDF input files are.
111  ! These can be overridden later in via the field manager?
112 
113  integer :: ntau, k,i,j,axes(3)
114  type(g_tracer_type), pointer :: g_tracer,g_tracer_next
115  character(len=fm_string_len) :: g_tracer_name,longname,units
116  real, dimension(:,:,:,:), pointer :: tr_field
117  real, dimension(:,:,:), pointer :: tr_ptr
118  real, dimension(HI%isd:HI%ied, HI%jsd:HI%jed,GV%ke) :: grid_tmask
119  integer, dimension(HI%isd:HI%ied, HI%jsd:HI%jed) :: grid_kmt
120 
121  register_mom_generic_tracer = .false.
122  if (associated(cs)) then
123  call mom_error(warning, "register_MOM_generic_tracer called with an "// &
124  "associated control structure.")
125  return
126  endif
127  allocate(cs)
128 
129 
130  !Register all the generic tracers used and create the list of them.
131  !This can be called by ALL PE's. No array fields allocated.
132  if (.not. g_registered) then
133  call generic_tracer_register()
134  g_registered = .true.
135  endif
136 
137 
138  ! Read all relevant parameters and write them to the model log.
139  call log_version(param_file, sub_name, version, "")
140  call get_param(param_file, sub_name, "GENERIC_TRACER_IC_FILE", cs%IC_file, &
141  "The file in which the generic trcer initial values can "//&
142  "be found, or an empty string for internal initialization.", &
143  default=" ")
144  if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
145  ! Add the directory if CS%IC_file is not already a complete path.
146  call get_param(param_file, sub_name, "INPUTDIR", inputdir, default=".")
147  cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
148  call log_param(param_file, sub_name, "INPUTDIR/GENERIC_TRACER_IC_FILE", cs%IC_file)
149  endif
150  call get_param(param_file, sub_name, "GENERIC_TRACER_IC_FILE_IS_Z", cs%Z_IC_file, &
151  "If true, GENERIC_TRACER_IC_FILE is in depth space, not "//&
152  "layer space.",default=.false.)
153  call get_param(param_file, sub_name, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
154  "If true, tracers may go through the initialization code "//&
155  "if they are not found in the restart files. Otherwise "//&
156  "it is a fatal error if tracers are not found in the "//&
157  "restart files of a restarted run.", default=.false.)
158 
159  cs%restart_CSp => restart_cs
160 
161 
162  ntau=1 ! MOM needs the fields at only one time step
163 
164 
165  ! At this point G%mask2dT and CS%diag%axesTL are not allocated.
166  ! postpone diag_registeration to initialize_MOM_generic_tracer
167 
168  !Fields cannot be diag registered as they are allocated and have to registered later.
169  grid_tmask(:,:,:) = 0.0
170  grid_kmt(:,:) = 0.0
171  axes(:) = -1
172 
173  !
174  ! Initialize all generic tracers
175  !
176  call generic_tracer_init(hi%isc,hi%iec,hi%jsc,hi%jec,hi%isd,hi%ied,hi%jsd,hi%jed,&
177  gv%ke,ntau,axes,grid_tmask,grid_kmt,set_time(0,0))
178 
179 
180  !
181  ! MOM-register the generic tracers
182  !
183 
184  !Get the tracer list
185  call generic_tracer_get_list(cs%g_tracer_list)
186  if (.NOT. associated(cs%g_tracer_list)) call mom_error(fatal, trim(sub_name)//&
187  ": No tracer in the list.")
188  ! For each tracer name get its T_prog index and get its fields
189 
190  g_tracer=>cs%g_tracer_list
191  do
192  call g_tracer_get_alias(g_tracer,g_tracer_name)
193 
194  call g_tracer_get_pointer(g_tracer,g_tracer_name,'field',tr_field)
195  call g_tracer_get_values(g_tracer,g_tracer_name,'longname', longname)
196  call g_tracer_get_values(g_tracer,g_tracer_name,'units',units )
197 
198  !!nnz: MOM field is 3D. Does this affect performance? Need it be override field?
199  tr_ptr => tr_field(:,:,:,1)
200  ! Register prognastic tracer for horizontal advection, diffusion, and restarts.
201  if (g_tracer_is_prog(g_tracer)) then
202  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
203  name=g_tracer_name, longname=longname, units=units, &
204  registry_diags=.false., & !### CHANGE TO TRUE?
205  restart_cs=restart_cs, mandatory=.not.cs%tracers_may_reinit)
206  else
207  call register_restart_field(tr_ptr, g_tracer_name, .not.cs%tracers_may_reinit, &
208  restart_cs, longname=longname, units=units)
209  endif
210 
211  !traverse the linked list till hit NULL
212  call g_tracer_get_next(g_tracer, g_tracer_next)
213  if (.NOT. associated(g_tracer_next)) exit
214  g_tracer=>g_tracer_next
215 
216  enddo
217 
218  register_mom_generic_tracer = .true.