MOM6
pseudo_salt_tracer Module Reference

Data Types

type  p3d
 
type  pseudo_salt_tracer_cs
 

Functions/Subroutines

logical function, public register_pseudo_salt_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 
subroutine, public initialize_pseudo_salt_tracer (restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp, tv)
 
subroutine, public pseudo_salt_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, tv, debug, evap_CFL_limit, minimum_forcing_depth)
 
integer function, public pseudo_salt_stock (h, stocks, G, GV, CS, names, units, stock_index)
 
subroutine, public pseudo_salt_tracer_surface_state (state, h, G, CS)
 
subroutine, public pseudo_salt_tracer_end (CS)
 

Variables

integer, parameter ntr_max = 1
 

Function/Subroutine Documentation

◆ initialize_pseudo_salt_tracer()

subroutine, public pseudo_salt_tracer::initialize_pseudo_salt_tracer ( logical, intent(in)  restart,
type(time_type), intent(in), target  day,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(pseudo_salt_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp,
type(diag_to_z_cs), pointer  diag_to_Z_CSp,
type(thermo_var_ptrs), intent(in)  tv 
)
Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)
[in]tvA structure pointing to various thermodynamic variables

Definition at line 203 of file pseudo_salt_tracer.F90.

References mom_tracer_registry::add_tracer_diagnostics(), mom_io::query_vardesc(), and mom_diag_to_z::register_z_tracer().

203  logical, intent(in) :: restart
204  type(time_type), target, intent(in) :: day
205  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
206  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
207  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
208  type(diag_ctrl), target, intent(in) :: diag
209  type(ocean_obc_type), pointer :: obc
210  type(pseudo_salt_tracer_cs), pointer :: cs
211  type(sponge_cs), pointer :: sponge_csp
212  type(diag_to_z_cs), pointer :: diag_to_z_csp
213  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
214 ! This subroutine initializes the CS%ntr tracer fields in tr(:,:,:,:)
215 ! and it sets up the tracer output.
216 
217 ! Arguments: restart - .true. if the fields have already been read from
218 ! a restart file.
219 ! (in) day - Time of the start of the run.
220 ! (in) G - The ocean's grid structure.
221 ! (in) GV - The ocean's vertical grid structure.
222 ! (in) h - Layer thickness, in m or kg m-2.
223 ! (in) diag - A structure that is used to regulate diagnostic output.
224 ! (in) OBC - This open boundary condition type specifies whether, where,
225 ! and what open boundary conditions are used.
226 ! (in/out) CS - The control structure returned by a previous call to
227 ! register_pseudo_salt_tracer.
228 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
229 ! they are in use. Otherwise this may be unassociated.
230 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
231 ! in depth space.
232  character(len=16) :: name ! A variable's name in a NetCDF file.
233  character(len=72) :: longname ! The long name of that variable.
234  character(len=48) :: units ! The dimensions of the variable.
235  character(len=48) :: flux_units ! The units for age tracer fluxes, either
236  ! years m3 s-1 or years kg s-1.
237  logical :: ok
238  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
239  integer :: isdb, iedb, jsdb, jedb
240 
241  if (.not.associated(cs)) return
242  if (cs%ntr < 1) return
243  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
244  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
245  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
246 
247  cs%Time => day
248  cs%diag => diag
249  name = "pseudo_salt"
250 
251  do m=1,cs%ntr
252  call query_vardesc(cs%tr_desc(m), name=name, caller="initialize_pseudo_salt_tracer")
253  if ((.not.restart) .or. (.not. &
254  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
255  do k=1,nz ; do j=jsd,jed ; do i=isd,ied
256  cs%tr(i,j,k,m) = tv%S(i,j,k)
257  enddo ; enddo ; enddo
258  endif
259  enddo ! Tracer loop
260 
261  if (associated(obc)) then
262  ! All tracers but the first have 0 concentration in their inflows. As this
263  ! is the default value, the following calls are unnecessary.
264  ! do m=1,CS%ntr
265  ! call add_tracer_OBC_values(trim(CS%tr_desc(m)%name), CS%tr_Reg, 0.0)
266  ! enddo
267  endif
268 
269  ! This needs to be changed if the units of tracer are changed above.
270  if (gv%Boussinesq) then ; flux_units = "g salt/(m^2 s)"
271  else ; flux_units = "g salt/(m^2 s)" ; endif
272 
273  do m=1,cs%ntr
274  ! Register the tracer for the restart file.
275  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
276  caller="initialize_pseudo_salt_tracer")
277  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
278  day, trim(longname) , trim(units))
279  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
280  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
281  trim(flux_units))
282  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
283  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
284  trim(flux_units))
285  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
286  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
287  trim(flux_units))
288  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
289  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
290  trim(flux_units))
291  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
292  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
293  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
294  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
295 
296 ! Register the tracer for horizontal advection & diffusion.
297  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
298  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
299  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
300  cs%tr_ady(m)%p,cs%tr_dfx(m)%p,cs%tr_dfy(m)%p)
301 
302  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
303  day, g, diag_to_z_csp)
304  enddo
305 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ pseudo_salt_stock()

integer function, public pseudo_salt_tracer::pseudo_salt_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(pseudo_salt_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units,
integer, intent(in), optional  stock_index 
)
Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)

Definition at line 422 of file pseudo_salt_tracer.F90.

References mom_io::query_vardesc().

Referenced by mom_tracer_flow_control::call_tracer_stocks().

422  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
423  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
424  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
425  real, dimension(:), intent(out) :: stocks
426  type(pseudo_salt_tracer_cs), pointer :: cs
427  character(len=*), dimension(:), intent(out) :: names
428  character(len=*), dimension(:), intent(out) :: units
429  integer, optional, intent(in) :: stock_index
430  integer :: pseudo_salt_stock
431 ! This function calculates the mass-weighted integral of all tracer stocks,
432 ! returning the number of stocks it has calculated. If the stock_index
433 ! is present, only the stock corresponding to that coded index is returned.
434 
435 ! Arguments: h - Layer thickness, in m or kg m-2.
436 ! (out) stocks - the mass-weighted integrated amount of each tracer,
437 ! in kg times concentration units.
438 ! (in) G - The ocean's grid structure.
439 ! (in) GV - The ocean's vertical grid structure.
440 ! (in) CS - The control structure returned by a previous call to
441 ! register_pseudo_salt_tracer.
442 ! (out) names - the names of the stocks calculated.
443 ! (out) units - the units of the stocks calculated.
444 ! (in,opt) stock_index - the coded index of a specific stock being sought.
445 ! Return value: the number of stocks calculated here.
446 
447  integer :: i, j, k, is, ie, js, je, nz, m
448  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
449 
450  pseudo_salt_stock = 0
451  if (.not.associated(cs)) return
452  if (cs%ntr < 1) return
453 
454  if (present(stock_index)) then ; if (stock_index > 0) then
455  ! Check whether this stock is available from this routine.
456 
457  ! No stocks from this routine are being checked yet. Return 0.
458  return
459  endif ; endif
460 
461  do m=1,1
462  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="pseudo_salt_stock")
463  units(m) = trim(units(m))//" kg"
464  stocks(m) = 0.0
465  do k=1,nz ; do j=js,je ; do i=is,ie
466  stocks(m) = stocks(m) + cs%diff(i,j,k,m) * &
467  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
468  enddo ; enddo ; enddo
469  stocks(m) = gv%H_to_kg_m2 * stocks(m)
470  enddo
471 
472  pseudo_salt_stock = cs%ntr
473 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pseudo_salt_tracer_column_physics()

subroutine, public pseudo_salt_tracer::pseudo_salt_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, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(pseudo_salt_tracer_cs), pointer  CS,
type(thermo_var_ptrs), intent(in)  tv,
logical, intent(in)  debug,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)
Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]dtThe amount of time covered by this call, in s
[in]tvA structure pointing to various thermodynamic variables

Definition at line 310 of file pseudo_salt_tracer.F90.

References mom_tracer_diabatic::applytracerboundaryfluxesinout(), and mom_tracer_diabatic::tracer_vertdiff().

310  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
311  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
312  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
313  type(forcing), intent(in) :: fluxes
314  real, intent(in) :: dt !< The amount of time covered by this call, in s
315  type(pseudo_salt_tracer_cs), pointer :: cs
316  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
317  logical, intent(in) :: debug
318  real, optional,intent(in) :: evap_cfl_limit
319  real, optional,intent(in) :: minimum_forcing_depth
320 
321 ! This subroutine applies diapycnal diffusion and any other column
322 ! tracer physics or chemistry to the tracers from this file.
323 ! This is a simple example of a set of advected passive tracers.
324 
325 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
326 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
327 ! (in) ea - an array to which the amount of fluid entrained
328 ! from the layer above during this call will be
329 ! added, in m or kg m-2.
330 ! (in) eb - an array to which the amount of fluid entrained
331 ! from the layer below during this call will be
332 ! added, in m or kg m-2.
333 ! (in) fluxes - A structure containing pointers to any possible
334 ! forcing fields. Unused fields have NULL ptrs.
335 ! (in) dt - The amount of time covered by this call, in s.
336 ! (in) G - The ocean's grid structure.
337 ! (in) GV - The ocean's vertical grid structure.
338 ! (in) CS - The control structure returned by a previous call to
339 ! register_pseudo_salt_tracer.
340 ! (in) tv - Thermodynamic structure with T and S
341 ! (in) evap_CFL_limit - Limits how much water can be fluxed out of the top layer
342 ! Stored previously in diabatic CS.
343 ! (in) minimum_forcing_depth - The smallest depth over which fluxes can be applied
344 ! Stored previously in diabatic CS.
345 ! (in) debug - Calculates checksums
346 !
347 ! The arguments to this subroutine are redundant in that
348 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
349 
350  real :: isecs_per_year = 1.0 / (365.0*86400.0)
351  real :: year, h_total, scale, htot, ih_limit
352  integer :: secs, days
353  integer :: i, j, k, is, ie, js, je, nz, m, k_max
354  real, allocatable :: local_tr(:,:,:)
355  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
356  real, dimension(:,:), pointer :: net_salt
357 
358  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
359  net_salt=>fluxes%netSalt
360 
361  if (.not.associated(cs)) return
362  if (cs%ntr < 1) return
363 
364  if (debug) then
365  call hchksum(tv%S,"salt pre pseudo-salt vertdiff", g%HI)
366  call hchksum(cs%tr(:,:,:,1),"pseudo_salt pre pseudo-salt vertdiff", g%HI)
367  endif
368 
369  ! This uses applyTracerBoundaryFluxesInOut, usually in ALE mode
370  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
371  do k=1,nz ;do j=js,je ; do i=is,ie
372  h_work(i,j,k) = h_old(i,j,k)
373  enddo ; enddo ; enddo;
374  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,1), dt, fluxes, h_work, &
375  evap_cfl_limit, minimum_forcing_depth, out_flux_optional=net_salt)
376  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
377  else
378  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
379  endif
380 
381  do k=1,nz ; do j=js,je ; do i=is,ie
382  cs%diff(i,j,k,1) = cs%tr(i,j,k,1)-tv%S(i,j,k)
383  enddo ; enddo ; enddo
384 
385  if(debug) then
386  call hchksum(tv%S,"salt post pseudo-salt vertdiff", g%HI)
387  call hchksum(cs%tr(:,:,:,1),"pseudo_salt post pseudo-salt vertdiff", g%HI)
388  endif
389 
390  allocate(local_tr(g%isd:g%ied,g%jsd:g%jed,nz))
391  do m=1,1
392  if (cs%id_tracer(m)>0) then
393  if (cs%mask_tracers) then
394  do k=1,nz ; do j=js,je ; do i=is,ie
395  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
396  local_tr(i,j,k) = cs%land_val(m)
397  else
398  local_tr(i,j,k) = cs%diff(i,j,k,m)
399  endif
400  enddo ; enddo ; enddo
401  else
402  do k=1,nz ; do j=js,je ; do i=is,ie
403  local_tr(i,j,k) = cs%tr(i,j,k,m)-tv%S(i,j,k)
404  enddo ; enddo ; enddo
405  endif ! CS%mask_tracers
406  call post_data(cs%id_tracer(m),local_tr,cs%diag)
407  endif ! CS%id_tracer(m)>0
408  if (cs%id_tr_adx(m)>0) &
409  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
410  if (cs%id_tr_ady(m)>0) &
411  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
412  if (cs%id_tr_dfx(m)>0) &
413  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
414  if (cs%id_tr_dfy(m)>0) &
415  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
416  enddo
417  deallocate(local_tr)
418 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ pseudo_salt_tracer_end()

subroutine, public pseudo_salt_tracer::pseudo_salt_tracer_end ( type(pseudo_salt_tracer_cs), pointer  CS)

Definition at line 506 of file pseudo_salt_tracer.F90.

Referenced by mom_tracer_flow_control::tracer_flow_control_end().

506  type(pseudo_salt_tracer_cs), pointer :: cs
507  integer :: m
508 
509  if (associated(cs)) then
510  if (associated(cs%tr)) deallocate(cs%tr)
511  if (associated(cs%tr)) deallocate(cs%diff)
512  do m=1,cs%ntr
513  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
514  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
515  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
516  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
517  enddo
518 
519  deallocate(cs)
520  endif
Here is the caller graph for this function:

◆ pseudo_salt_tracer_surface_state()

subroutine, public pseudo_salt_tracer::pseudo_salt_tracer_surface_state ( type(surface), intent(inout)  state,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(pseudo_salt_tracer_cs), pointer  CS 
)
Parameters
[in]gThe ocean's grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)

Definition at line 477 of file pseudo_salt_tracer.F90.

References coupler_util::set_coupler_values().

477  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
478  type(surface), intent(inout) :: state
479  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
480  type(pseudo_salt_tracer_cs), pointer :: cs
481 ! This particular tracer package does not report anything back to the coupler.
482 ! The code that is here is just a rough guide for packages that would.
483 ! Arguments: state - A structure containing fields that describe the
484 ! surface state of the ocean.
485 ! (in) h - Layer thickness, in m or kg m-2.
486 ! (in) G - The ocean's grid structure.
487 ! (in) CS - The control structure returned by a previous call to
488 ! register_pseudo_salt_tracer.
489  integer :: m, is, ie, js, je
490  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
491 
492  if (.not.associated(cs)) return
493 
494  if (cs%coupled_tracers) then
495  do m=1,cs%ntr
496  ! This call loads the surface vlues into the appropriate array in the
497  ! coupler-type structure.
498  call set_coupler_values(cs%tr(:,:,1,m), state%tr_fields, cs%ind_tr(m), &
499  ind_csurf, is, ie, js, je)
500  enddo
501  endif
502 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:

◆ register_pseudo_salt_tracer()

logical function, public pseudo_salt_tracer::register_pseudo_salt_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(pseudo_salt_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)
Parameters
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters

Definition at line 130 of file pseudo_salt_tracer.F90.

References atmos_ocean_fluxes_mod::aof_set_coupler_flux(), mom_error_handler::mom_error(), ntr_max, mom_io::query_vardesc(), and mom_io::var_desc().

130  type(hor_index_type), intent(in) :: hi
131  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
132  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
133  type(pseudo_salt_tracer_cs), pointer :: cs
134  type(tracer_registry_type), pointer :: tr_reg
135  type(mom_restart_cs), pointer :: restart_cs
136 ! This subroutine is used to register tracer fields and subroutines
137 ! to be used with MOM.
138 ! Arguments: HI - A horizontal index type structure.
139 ! (in) GV - The ocean's vertical grid structure.
140 ! (in) param_file - A structure indicating the open file to parse for
141 ! model parameter values.
142 ! (in/out) CS - A pointer that is set to point to the control structure
143 ! for this module
144 ! (in/out) tr_Reg - A pointer that is set to point to the control structure
145 ! for the tracer advection and diffusion module.
146 ! (in) restart_CS - A pointer to the restart control structure.
147 
148 ! This include declares and sets the variable "version".
149 #include "version_variable.h"
150  character(len=40) :: mdl = "pseudo_salt_tracer" ! This module's name.
151  character(len=200) :: inputdir ! The directory where the input files are.
152  character(len=48) :: var_name ! The variable's name.
153  character(len=3) :: name_tag ! String for creating identifying pseudo_salt
154  real, pointer :: tr_ptr(:,:,:) => null()
155  logical :: register_pseudo_salt_tracer
156  integer :: isd, ied, jsd, jed, nz, m, i, j
157  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
158 
159  if (associated(cs)) then
160  call mom_error(warning, "register_pseudo_salt_tracer called with an "// &
161  "associated control structure.")
162  return
163  endif
164  allocate(cs)
165 
166  ! Read all relevant parameters and write them to the model log.
167  call log_version(param_file, mdl, version, "")
168 
169  cs%ntr = ntr_max
170  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
171  allocate(cs%diff(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%diff(:,:,:,:) = 0.0
172 
173  do m=1,cs%ntr
174  ! This is needed to force the compiler not to do a copy in the registration
175  ! calls. Curses on the designers and implementers of Fortran90.
176  cs%tr_desc(m) = var_desc(trim("pseudo_salt_diff"), "kg", &
177  "Difference between pseudo salt passive tracer and salt tracer", caller=mdl)
178  tr_ptr => cs%tr(:,:,:,m)
179  call query_vardesc(cs%tr_desc(m), name=var_name, caller="register_pseudo_salt_tracer")
180  ! Register the tracer for the restart file.
181  call register_restart_field(tr_ptr, cs%tr_desc(m), &
182  .not. cs%pseudo_salt_may_reinit, restart_cs)
183  ! Register the tracer for horizontal advection & diffusion.
184  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
185  tr_desc_ptr=cs%tr_desc(m))
186 
187  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
188  ! values to the coupler (if any). This is meta-code and its arguments will
189  ! currently (deliberately) give fatal errors if it is used.
190  if (cs%coupled_tracers) &
191  cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//'_flux', &
192  flux_type=' ', implementation=' ', caller="register_pseudo_salt_tracer")
193  enddo
194 
195  cs%tr_Reg => tr_reg
196  cs%restart_CSp => restart_cs
197  register_pseudo_salt_tracer = .true.
198 
Here is the call graph for this function:

Variable Documentation

◆ ntr_max

integer, parameter pseudo_salt_tracer::ntr_max = 1
private

Definition at line 88 of file pseudo_salt_tracer.F90.

Referenced by register_pseudo_salt_tracer().

88 integer, parameter :: ntr_max = 1