MOM6
pseudo_salt_tracer Module Reference

Detailed Description

A tracer package that mimics salinity.

By Andrew Shao, 2016

This file contains the routines necessary to model a passive tracer that uses the same boundary fluxes as salinity. At the beginning of the run, salt is set to the same as tvS. Any deviations between this salt-like tracer and tvS signifies a difference between how active and passive tracers are treated.

Data Types

type  pseudo_salt_tracer_cs
 The control structure for the pseudo-salt tracer. More...
 

Functions/Subroutines

logical function, public register_pseudo_salt_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Register the pseudo-salt tracer with MOM6. More...
 
subroutine, public initialize_pseudo_salt_tracer (restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, tv)
 Initialize the pseudo-salt tracer. More...
 
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)
 Apply sources, sinks and diapycnal diffusion to the tracers in this package. More...
 
integer function, public pseudo_salt_stock (h, stocks, G, GV, CS, names, units, stock_index)
 Calculates the mass-weighted integral of all tracer stocks, 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...
 
subroutine, public pseudo_salt_tracer_surface_state (state, h, G, CS)
 This subroutine extracts the surface fields from this tracer package that are to be shared with the atmosphere in coupled configurations. This particular tracer package does not report anything back to the coupler. More...
 
subroutine, public pseudo_salt_tracer_end (CS)
 Deallocate memory associated with this tracer package. More...
 

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( g %isd: g %ied, g %jsd: g %jed, g %ke), 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(thermo_var_ptrs), intent(in)  tv 
)

Initialize the pseudo-salt tracer.

Parameters
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]diagA structure that is used to regulate diagnostic output.
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used.
csThe control structure returned by a previous call to register_pseudo_salt_tracer.
sponge_cspPointer to the control structure for the sponges.
[in]tvA structure pointing to various thermodynamic variables

Definition at line 117 of file pseudo_salt_tracer.F90.

117  logical, intent(in) :: restart !< .true. if the fields have already
118  !! been read from a restart file.
119  type(time_type), target, intent(in) :: day !< Time of the start of the run.
120  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
121  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
122  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
123  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
124  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
125  !! diagnostic output.
126  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
127  !! whether, where, and what open boundary
128  !! conditions are used.
129  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
130  !! call to register_pseudo_salt_tracer.
131  type(sponge_CS), pointer :: sponge_CSp !< Pointer to the control structure for the sponges.
132  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
133 ! This subroutine initializes the tracer fields in CS%ps(:,:,:).
134 
135  ! Local variables
136  character(len=16) :: name ! A variable's name in a NetCDF file.
137  character(len=72) :: longname ! The long name of that variable.
138  character(len=48) :: units ! The dimensions of the variable.
139  character(len=48) :: flux_units ! The units for age tracer fluxes, either
140  ! years m3 s-1 or years kg s-1.
141  logical :: OK
142  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
143  integer :: IsdB, IedB, JsdB, JedB
144 
145  if (.not.associated(cs)) return
146  if (.not.associated(cs%diff)) return
147 
148  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
149  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
150  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
151 
152  cs%Time => day
153  cs%diag => diag
154  name = "pseudo_salt"
155 
156  call query_vardesc(cs%tr_desc, name=name, caller="initialize_pseudo_salt_tracer")
157  if ((.not.restart) .or. (.not.query_initialized(cs%ps, name, cs%restart_CSp))) then
158  do k=1,nz ; do j=jsd,jed ; do i=isd,ied
159  cs%ps(i,j,k) = tv%S(i,j,k)
160  enddo ; enddo ; enddo
161  endif
162 
163  if (associated(obc)) then
164  ! Steal from updated DOME in the fullness of time.
165  endif
166 
167  cs%id_psd = register_diag_field("ocean_model", "pseudo_salt_diff", cs%diag%axesTL, &
168  day, "Difference between pseudo salt passive tracer and salt tracer", "psu")
169 

◆ pseudo_salt_stock()

integer function, public pseudo_salt_tracer::pseudo_salt_stock ( real, dimension(szi_(g),szj_(g),szk_(g)), 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 
)

Calculates the mass-weighted integral of all tracer stocks, 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].
csThe control structure returned by a previous call to register_pseudo_salt_tracer.
[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 252 of file pseudo_salt_tracer.F90.

252  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
253  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
254  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
255  real, dimension(:), intent(out) :: stocks !< the mass-weighted integrated amount of each
256  !! tracer, in kg times concentration units [kg conc].
257  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
258  !! call to register_pseudo_salt_tracer.
259  character(len=*), dimension(:), intent(out) :: names !< The names of the stocks calculated.
260  character(len=*), dimension(:), intent(out) :: units !< The units of the stocks calculated.
261  integer, optional, intent(in) :: stock_index !< The coded index of a specific stock
262  !! being sought.
263  integer :: pseudo_salt_stock !< Return value: the number of
264  !! stocks calculated here.
265 
266 ! This function calculates the mass-weighted integral of all tracer stocks,
267 ! returning the number of stocks it has calculated. If the stock_index
268 ! is present, only the stock corresponding to that coded index is returned.
269 
270  integer :: i, j, k, is, ie, js, je, nz
271  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
272 
273  pseudo_salt_stock = 0
274  if (.not.associated(cs)) return
275  if (.not.associated(cs%diff)) return
276 
277  if (present(stock_index)) then ; if (stock_index > 0) then
278  ! Check whether this stock is available from this routine.
279 
280  ! No stocks from this routine are being checked yet. Return 0.
281  return
282  endif ; endif
283 
284  call query_vardesc(cs%tr_desc, name=names(1), units=units(1), caller="pseudo_salt_stock")
285  units(1) = trim(units(1))//" kg"
286  stocks(1) = 0.0
287  do k=1,nz ; do j=js,je ; do i=is,ie
288  stocks(1) = stocks(1) + cs%diff(i,j,k) * &
289  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
290  enddo ; enddo ; enddo
291  stocks(1) = gv%H_to_kg_m2 * stocks(1)
292 
293  pseudo_salt_stock = 1
294 

◆ pseudo_salt_tracer_column_physics()

subroutine, public pseudo_salt_tracer::pseudo_salt_tracer_column_physics ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_old,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_new,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  ea,
real, dimension(szi_(g),szj_(g),szk_(g)), 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 
)

Apply sources, sinks and diapycnal diffusion to the tracers in this package.

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]eaan array to which the amount of fluid entrained
[in]eban array to which the amount of fluid entrained
[in]fluxesA structure containing pointers to thermodynamic and tracer forcing fields. Unused fields have NULL ptrs.
[in]dtThe amount of time covered by this call [s]
csThe control structure returned by a previous call to register_pseudo_salt_tracer.
[in]tvA structure pointing to various thermodynamic variables
[in]debugIf true calculate checksums
[in]evap_cfl_limitLimit on the fraction of the water that can be fluxed out of the top layer in a timestep [nondim]
[in]minimum_forcing_depthThe smallest depth over which fluxes can be applied [m]

Definition at line 175 of file pseudo_salt_tracer.F90.

175  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
176  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
177  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
178  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
179  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
180  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
181  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
182  intent(in) :: ea !< an array to which the amount of fluid entrained
183  !! from the layer above during this call will be
184  !! added [H ~> m or kg m-2].
185  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
186  intent(in) :: eb !< an array to which the amount of fluid entrained
187  !! from the layer below during this call will be
188  !! added [H ~> m or kg m-2].
189  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
190  !! and tracer forcing fields. Unused fields have NULL ptrs.
191  real, intent(in) :: dt !< The amount of time covered by this call [s]
192  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
193  !! call to register_pseudo_salt_tracer.
194  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
195  logical, intent(in) :: debug !< If true calculate checksums
196  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
197  !! be fluxed out of the top layer in a timestep [nondim]
198  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
199  !! fluxes can be applied [m]
200 
201 ! This subroutine applies diapycnal diffusion and any other column
202 ! tracer physics or chemistry to the tracers from this file.
203 
204 ! The arguments to this subroutine are redundant in that
205 ! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
206 
207  ! Local variables
208  real :: year, h_total, scale, htot, Ih_limit
209  integer :: secs, days
210  integer :: i, j, k, is, ie, js, je, nz, k_max
211  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
212 
213  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
214 
215  if (.not.associated(cs)) return
216  if (.not.associated(cs%diff)) return
217 
218  if (debug) then
219  call hchksum(tv%S,"salt pre pseudo-salt vertdiff", g%HI)
220  call hchksum(cs%ps,"pseudo_salt pre pseudo-salt vertdiff", g%HI)
221  endif
222 
223  ! This uses applyTracerBoundaryFluxesInOut, usually in ALE mode
224  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
225  do k=1,nz ; do j=js,je ; do i=is,ie
226  h_work(i,j,k) = h_old(i,j,k)
227  enddo ; enddo ; enddo
228  call applytracerboundaryfluxesinout(g, gv, cs%ps, dt, fluxes, h_work, &
229  evap_cfl_limit, minimum_forcing_depth, out_flux_optional=fluxes%netSalt)
230  call tracer_vertdiff(h_work, ea, eb, dt, cs%ps, g, gv)
231  else
232  call tracer_vertdiff(h_old, ea, eb, dt, cs%ps, g, gv)
233  endif
234 
235  do k=1,nz ; do j=js,je ; do i=is,ie
236  cs%diff(i,j,k) = cs%ps(i,j,k)-tv%S(i,j,k)
237  enddo ; enddo ; enddo
238 
239  if (debug) then
240  call hchksum(tv%S,"salt post pseudo-salt vertdiff", g%HI)
241  call hchksum(cs%ps,"pseudo_salt post pseudo-salt vertdiff", g%HI)
242  endif
243 
244  if (cs%id_psd>0) call post_data(cs%id_psd, cs%diff, cs%diag)
245 

◆ pseudo_salt_tracer_end()

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

Deallocate memory associated with this tracer package.

Parameters
csThe control structure returned by a previous call to register_pseudo_salt_tracer.

Definition at line 324 of file pseudo_salt_tracer.F90.

324  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
325  !! call to register_pseudo_salt_tracer.
326  integer :: m
327 
328  if (associated(cs)) then
329  if (associated(cs%ps)) deallocate(cs%ps)
330  if (associated(cs%diff)) deallocate(cs%diff)
331  deallocate(cs)
332  endif

◆ pseudo_salt_tracer_surface_state()

subroutine, public pseudo_salt_tracer::pseudo_salt_tracer_surface_state ( type(surface), intent(inout)  state,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(pseudo_salt_tracer_cs), pointer  CS 
)

This subroutine extracts the surface fields from this tracer package that are to be shared with the atmosphere in coupled configurations. This particular tracer package does not report anything back to the coupler.

Parameters
[in]gThe ocean's grid structure.
[in,out]stateA structure containing fields that describe the surface state of the ocean.
[in]hLayer thickness [H ~> m or kg m-2].
csThe control structure returned by a previous call to register_pseudo_salt_tracer.

Definition at line 301 of file pseudo_salt_tracer.F90.

301  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
302  type(surface), intent(inout) :: state !< A structure containing fields that
303  !! describe the surface state of the ocean.
304  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
305  intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
306  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
307  !! call to register_pseudo_salt_tracer.
308 
309  ! This particular tracer package does not report anything back to the coupler.
310  ! The code that is here is just a rough guide for packages that would.
311 
312  integer :: m, is, ie, js, je, isd, ied, jsd, jed
313  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
314  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
315 
316  if (.not.associated(cs)) return
317 
318  ! By design, this tracer package does not return any surface states.
319 

◆ 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 
)

Register the pseudo-salt tracer with MOM6.

Parameters
[in]hiA horizontal index type structure
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters
csThe control structure returned by a previous call to register_pseudo_salt_tracer.
tr_regA pointer that is set to point to the control structure for the tracer advection and diffusion module
restart_csA pointer to the restart control structure

Definition at line 60 of file pseudo_salt_tracer.F90.

60  type(hor_index_type), intent(in) :: HI !< A horizontal index type structure
61  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
62  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
63  type(pseudo_salt_tracer_CS), pointer :: CS !< The control structure returned by a previous
64  !! call to register_pseudo_salt_tracer.
65  type(tracer_registry_type), pointer :: tr_Reg !< A pointer that is set to point to the control
66  !! structure for the tracer advection and
67  !! diffusion module
68  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure
69 ! This subroutine is used to register tracer fields and subroutines
70 ! to be used with MOM.
71 
72  ! Local variables
73  character(len=40) :: mdl = "pseudo_salt_tracer" ! This module's name.
74  character(len=200) :: inputdir ! The directory where the input files are.
75  character(len=48) :: var_name ! The variable's name.
76  character(len=3) :: name_tag ! String for creating identifying pseudo_salt
77 ! This include declares and sets the variable "version".
78 #include "version_variable.h"
79  real, pointer :: tr_ptr(:,:,:) => null()
80  logical :: register_pseudo_salt_tracer
81  integer :: isd, ied, jsd, jed, nz, i, j
82  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
83 
84  if (associated(cs)) then
85  call mom_error(warning, "register_pseudo_salt_tracer called with an "// &
86  "associated control structure.")
87  return
88  endif
89  allocate(cs)
90 
91  ! Read all relevant parameters and write them to the model log.
92  call log_version(param_file, mdl, version, "")
93 
94  allocate(cs%ps(isd:ied,jsd:jed,nz)) ; cs%ps(:,:,:) = 0.0
95  allocate(cs%diff(isd:ied,jsd:jed,nz)) ; cs%diff(:,:,:) = 0.0
96 
97  cs%tr_desc = var_desc(trim("pseudo_salt"), "psu", &
98  "Pseudo salt passive tracer", caller=mdl)
99 
100  tr_ptr => cs%ps(:,:,:)
101  call query_vardesc(cs%tr_desc, name=var_name, caller="register_pseudo_salt_tracer")
102  ! Register the tracer for horizontal advection, diffusion, and restarts.
103  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, name="pseudo_salt", &
104  longname="Pseudo salt passive tracer", units="psu", &
105  registry_diags=.true., restart_cs=restart_cs, &
106  mandatory=.not.cs%pseudo_salt_may_reinit)
107 
108  cs%tr_Reg => tr_reg
109  cs%restart_CSp => restart_cs
110  register_pseudo_salt_tracer = .true.
111