MOM6
mom_tracer_z_init Module Reference

Detailed Description

Used to initialize tracers from a depth- (or z*-) space file.

Functions/Subroutines

logical function, public tracer_z_init (tr, h, filename, tr_name, G, US, missing_val, land_val)
 This function initializes a tracer by reading a Z-space file, returning .true. if this appears to have been successful, and false otherwise. More...
 
subroutine, public tracer_z_init_array (tr_in, z_edges, nk_data, e, land_fill, G, nlay, nlevs, eps_z, tr)
 Layer model routine for remapping tracers from pseudo-z coordinates into layers defined by target interface positions. More...
 
subroutine read_z_edges (filename, tr_name, z_edges, nz_out, has_edges, use_missing, missing, scale)
 This subroutine reads the vertical coordinate data for a field from a NetCDF file. It also might read the missing value attribute for that same field. More...
 
subroutine find_overlap (e, Z_top, Z_bot, k_max, k_start, k_top, k_bot, wt, z1, z2)
 Determines the layers bounded by interfaces e that overlap with the depth range between Z_top and Z_bot, and the fractional weights of each layer. It also calculates the normalized relative depths of the range of each layer that overlaps that depth range. More...
 
real function find_limited_slope (val, e, k)
 This subroutine determines a limited slope for val to be advected with a piecewise limited scheme. More...
 
subroutine, public determine_temperature (temp, salt, R_tgt, p_ref, niter, land_fill, h, k_start, G, US, eos, h_massless)
 This subroutine determines the potential temperature and salinity that is consistent with the target density using provided initial guess. More...
 

Function/Subroutine Documentation

◆ determine_temperature()

subroutine, public mom_tracer_z_init::determine_temperature ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  temp,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  salt,
real, dimension( g %ke), intent(in)  R_tgt,
real, intent(in)  p_ref,
integer, intent(in)  niter,
real, intent(in)  land_fill,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
integer, intent(in)  k_start,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(eos_type), pointer  eos,
real, intent(in), optional  h_massless 
)

This subroutine determines the potential temperature and salinity that is consistent with the target density using provided initial guess.

Parameters
[in]gThe ocean's grid structure
[in,out]temppotential temperature [degC]
[in,out]saltsalinity [PSU]
[in]r_tgtdesired potential density [R ~> kg m-3].
[in]p_refreference pressure [R L2 T-2 ~> Pa].
[in]nitermaximum number of iterations
[in]k_startstarting index (i.e. below the buffer layer)
[in]land_fillland fill value
[in]hlayer thickness, used only to avoid working on
[in]usA dimensional unit scaling type
eosseawater equation of state control structure
[in]h_masslessA threshold below which a layer is determined to be massless [H ~> m or kg m-2]

Definition at line 614 of file MOM_tracer_Z_init.F90.

614  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
615  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
616  intent(inout) :: temp !< potential temperature [degC]
617  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
618  intent(inout) :: salt !< salinity [PSU]
619  real, dimension(SZK_(G)), intent(in) :: r_tgt !< desired potential density [R ~> kg m-3].
620  real, intent(in) :: p_ref !< reference pressure [R L2 T-2 ~> Pa].
621  integer, intent(in) :: niter !< maximum number of iterations
622  integer, intent(in) :: k_start !< starting index (i.e. below the buffer layer)
623  real, intent(in) :: land_fill !< land fill value
624  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
625  intent(in) :: h !< layer thickness, used only to avoid working on
626  !! massless layers [H ~> m or kg m-2]
627  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
628  type(eos_type), pointer :: eos !< seawater equation of state control structure
629  real, optional, intent(in) :: h_massless !< A threshold below which a layer is
630  !! determined to be massless [H ~> m or kg m-2]
631 
632  real, parameter :: t_max = 31.0, t_min = -2.0
633  ! Local variables (All of which need documentation!)
634  real, dimension(SZI_(G),SZK_(G)) :: &
635  t, s, dt, ds, &
636  rho, & ! Layer densities [R ~> kg m-3]
637  hin, & ! Input layer thicknesses [H ~> m or kg m-2]
638  drho_dt, & ! Partial derivative of density with temperature [R degC-1 ~> kg m-3 degC-1]
639  drho_ds ! Partial derivative of density with salinity [R ppt-1 ~> kg m-3 ppt-1]
640  real, dimension(SZI_(G)) :: press ! Reference pressures [R L2 T-2 ~> Pa]
641  real :: dt_ds_gauge ! The relative penalizing of temperature to salinity changes when
642  ! minimizing property changes while correcting density [degC ppt-1].
643  real :: i_denom ! The inverse of the magnitude squared of the density gradient in
644  ! T-S space streched with dT_dS_gauge [ppt2 R-2 ~> ppt2 m6 kg-2]
645  logical :: adjust_salt, old_fit
646  real :: s_min, s_max
647  real :: tol_t ! The tolerance for temperature matches [degC]
648  real :: tol_s ! The tolerance for salinity matches [ppt]
649  real :: tol_rho ! The tolerance for density matches [R ~> kg m-3]
650  real :: max_t_adj, max_s_adj
651  integer, dimension(2) :: eosdom ! The i-computational domain for the equation of state
652  integer :: i, j, k, kz, is, ie, js, je, nz, itt
653 
654  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
655 
656  ! These hard coded parameters need to be set properly.
657  s_min = 0.5 ; s_max = 65.0
658  max_t_adj = 1.0 ; max_s_adj = 0.5
659  tol_t=1.e-4 ; tol_s=1.e-4 ; tol_rho = 1.e-4*us%kg_m3_to_R
660  old_fit = .true. ! reproduces siena behavior
661 
662  ! ### The whole determine_temperature subroutine needs to be reexamined, both the algorithms
663  ! and the extensive use of hard-coded dimensional parameters.
664 
665  ! We will switch to the newer method which simultaneously adjusts
666  ! temp and salt based on the ratio of the thermal and haline coefficients, once it is tested.
667 
668  press(:) = p_ref
669  eosdom(:) = eos_domain(g%HI)
670 
671  do j=js,je
672  ds(:,:) = 0. ! Needs to be zero everywhere since there is a maxval(abs(dS)) later...
673  t(:,:) = temp(:,j,:)
674  s(:,:) = salt(:,j,:)
675  hin(:,:) = h(:,j,:)
676  dt(:,:) = 0.0
677  adjust_salt = .true.
678  iter_loop: do itt = 1,niter
679  do k=1,nz
680  call calculate_density(t(:,k), s(:,k), press, rho(:,k), eos, eosdom )
681  call calculate_density_derivs(t(:,k), s(:,k), press, drho_dt(:,k), drho_ds(:,k), &
682  eos, eosdom )
683  enddo
684  do k=k_start,nz ; do i=is,ie
685 ! if (abs(rho(i,k)-R_tgt(k))>tol_rho .and. hin(i,k)>h_massless .and. abs(T(i,k)-land_fill) < epsln) then
686  if (abs(rho(i,k)-r_tgt(k))>tol_rho) then
687  if (old_fit) then
688  dt(i,k) = max(min((r_tgt(k)-rho(i,k)) / drho_dt(i,k), max_t_adj), -max_t_adj)
689  t(i,k) = max(min(t(i,k)+dt(i,k), t_max), t_min)
690  else
691  dt_ds_gauge = 10.0 ! 10 degC is weighted equivalently to 1 ppt.
692  i_denom = 1.0 / (drho_ds(i,k)**2 + dt_ds_gauge**2*drho_dt(i,k)**2)
693  ds(i,k) = (r_tgt(k)-rho(i,k)) * drho_ds(i,k) * i_denom
694  dt(i,k) = (r_tgt(k)-rho(i,k)) * dt_ds_gauge**2*drho_dt(i,k) * i_denom
695 
696  t(i,k) = max(min(t(i,k)+dt(i,k), t_max), t_min)
697  s(i,k) = max(min(s(i,k)+ds(i,k), s_max), s_min)
698  endif
699  endif
700  enddo ; enddo
701  if (maxval(abs(dt)) < tol_t) then
702  adjust_salt = .false.
703  exit iter_loop
704  endif
705  enddo iter_loop
706 
707  if (adjust_salt .and. old_fit) then ; do itt = 1,niter
708  do k=1,nz
709  call calculate_density(t(:,k), s(:,k), press, rho(:,k), eos, eosdom )
710  call calculate_density_derivs(t(:,k), s(:,k), press, drho_dt(:,k), drho_ds(:,k), &
711  eos, eosdom )
712  enddo
713  do k=k_start,nz ; do i=is,ie
714 ! if (abs(rho(i,k)-R_tgt(k))>tol_rho .and. hin(i,k)>h_massless .and. abs(T(i,k)-land_fill) < epsln ) then
715  if (abs(rho(i,k)-r_tgt(k)) > tol_rho) then
716  ds(i,k) = max(min((r_tgt(k)-rho(i,k)) / drho_ds(i,k), max_s_adj), -max_s_adj)
717  s(i,k) = max(min(s(i,k)+ds(i,k), s_max), s_min)
718  endif
719  enddo ; enddo
720  if (maxval(abs(ds)) < tol_s) exit
721  enddo ; endif
722 
723  temp(:,j,:) = t(:,:)
724  salt(:,j,:) = s(:,:)
725  enddo
726 

◆ find_limited_slope()

real function mom_tracer_z_init::find_limited_slope ( real, dimension(:), intent(in)  val,
real, dimension(:), intent(in)  e,
integer, intent(in)  k 
)
private

This subroutine determines a limited slope for val to be advected with a piecewise limited scheme.

Parameters
[in]valAn column the values that are being interpolated.
[in]eA column's interface heights [Z ~> m] or other units.
[in]kThe layer whose slope is being determined.
Returns
The normalized slope in the intracell distribution of val.

Definition at line 580 of file MOM_tracer_Z_init.F90.

580  real, dimension(:), intent(in) :: val !< An column the values that are being interpolated.
581  real, dimension(:), intent(in) :: e !< A column's interface heights [Z ~> m] or other units.
582  integer, intent(in) :: k !< The layer whose slope is being determined.
583  real :: slope !< The normalized slope in the intracell distribution of val.
584  ! Local variables
585  real :: amn, cmn
586  real :: d1, d2
587 
588  if ((val(k)-val(k-1)) * (val(k)-val(k+1)) >= 0.0) then
589  slope = 0.0 ! ; curvature = 0.0
590  else
591  d1 = 0.5*(e(k-1)-e(k+1)) ; d2 = 0.5*(e(k)-e(k+2))
592  if (d1*d2 > 0.0) then
593  slope = ((d1**2)*(val(k+1) - val(k)) + (d2**2)*(val(k) - val(k-1))) * &
594  (e(k) - e(k+1)) / (d1*d2*(d1+d2))
595  ! slope = 0.5*(val(k+1) - val(k-1))
596  ! This is S.J. Lin's form of the PLM limiter.
597  amn = min(abs(slope), 2.0*(max(val(k-1), val(k), val(k+1)) - val(k)))
598  cmn = 2.0*(val(k) - min(val(k-1), val(k), val(k+1)))
599  slope = sign(1.0, slope) * min(amn, cmn)
600 
601  ! min(abs(slope), 2.0*(max(val(k-1),val(k),val(k+1)) - val(k)), &
602  ! 2.0*(val(k) - min(val(k-1),val(k),val(k+1))))
603  ! curvature = 0.0
604  else
605  slope = 0.0 ! ; curvature = 0.0
606  endif
607  endif
608 

◆ find_overlap()

subroutine mom_tracer_z_init::find_overlap ( real, dimension(:), intent(in)  e,
real, intent(in)  Z_top,
real, intent(in)  Z_bot,
integer, intent(in)  k_max,
integer, intent(in)  k_start,
integer, intent(out)  k_top,
integer, intent(out)  k_bot,
real, dimension(:), intent(out)  wt,
real, dimension(:), intent(out)  z1,
real, dimension(:), intent(out)  z2 
)
private

Determines the layers bounded by interfaces e that overlap with the depth range between Z_top and Z_bot, and the fractional weights of each layer. It also calculates the normalized relative depths of the range of each layer that overlaps that depth range.

Parameters
[in]eColumn interface heights, [Z ~> m] or other units.
[in]z_topTop of range being mapped to, in the units of e [Z ~> m].
[in]z_botBottom of range being mapped to, in the units of e [Z ~> m].
[in]k_maxNumber of valid layers.
[in]k_startLayer at which to start searching.
[out]k_topIndices of top layers that overlap with the depth range.
[out]k_botIndices of bottom layers that overlap with the depth range.
[out]wtRelative weights of each layer from k_top to k_bot [nondim].
[out]z1Depth of the top limits of the part of a layer that contributes to a depth level, relative to the cell center and normalized by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.
[out]z2Depths of the bottom limit of the part of a layer that contributes to a depth level, relative to the cell center and normalized by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.

Definition at line 518 of file MOM_tracer_Z_init.F90.

518  real, dimension(:), intent(in) :: e !< Column interface heights, [Z ~> m] or other units.
519  real, intent(in) :: z_top !< Top of range being mapped to, in the units of e [Z ~> m].
520  real, intent(in) :: z_bot !< Bottom of range being mapped to, in the units of e [Z ~> m].
521  integer, intent(in) :: k_max !< Number of valid layers.
522  integer, intent(in) :: k_start !< Layer at which to start searching.
523  integer, intent(out) :: k_top !< Indices of top layers that overlap with the depth range.
524  integer, intent(out) :: k_bot !< Indices of bottom layers that overlap with the depth range.
525  real, dimension(:), intent(out) :: wt !< Relative weights of each layer from k_top to k_bot [nondim].
526  real, dimension(:), intent(out) :: z1 !< Depth of the top limits of the part of
527  !! a layer that contributes to a depth level, relative to the cell center and normalized
528  !! by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.
529  real, dimension(:), intent(out) :: z2 !< Depths of the bottom limit of the part of
530  !! a layer that contributes to a depth level, relative to the cell center and normalized
531  !! by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.
532  ! Local variables
533  real :: ih, e_c, tot_wt, i_totwt
534  integer :: k
535 
536  wt(:) = 0.0 ; z1(:) = 0.0 ; z2(:) = 0.0 ; k_bot = k_max
537 
538  do k=k_start,k_max ; if (e(k+1) < z_top) exit ; enddo
539  k_top = k
540  if (k_top > k_max) return
541 
542  ! Determine the fractional weights of each layer.
543  ! Note that by convention, e and Z_int decrease with increasing k.
544  if (e(k+1) <= z_bot) then
545  wt(k) = 1.0 ; k_bot = k
546  ih = 0.0 ; if (e(k) /= e(k+1)) ih = 1.0 / (e(k)-e(k+1))
547  e_c = 0.5*(e(k)+e(k+1))
548  z1(k) = (e_c - min(e(k), z_top)) * ih
549  z2(k) = (e_c - z_bot) * ih
550  else
551  wt(k) = min(e(k),z_top) - e(k+1) ; tot_wt = wt(k) ! These are always > 0.
552  if (e(k) /= e(k+1)) then
553  z1(k) = (0.5*(e(k)+e(k+1)) - min(e(k), z_top)) / (e(k)-e(k+1))
554  else ; z1(k) = -0.5 ; endif
555  z2(k) = 0.5
556  k_bot = k_max
557  do k=k_top+1,k_max
558  if (e(k+1) <= z_bot) then
559  k_bot = k
560  wt(k) = e(k) - z_bot ; z1(k) = -0.5
561  if (e(k) /= e(k+1)) then
562  z2(k) = (0.5*(e(k)+e(k+1)) - z_bot) / (e(k)-e(k+1))
563  else ; z2(k) = 0.5 ; endif
564  else
565  wt(k) = e(k) - e(k+1) ; z1(k) = -0.5 ; z2(k) = 0.5
566  endif
567  tot_wt = tot_wt + wt(k) ! wt(k) is always > 0.
568  if (k>=k_bot) exit
569  enddo
570 
571  i_totwt = 0.0 ; if (tot_wt > 0.0) i_totwt = 1.0 / tot_wt
572  do k=k_top,k_bot ; wt(k) = i_totwt*wt(k) ; enddo
573  endif
574 

◆ read_z_edges()

subroutine mom_tracer_z_init::read_z_edges ( character(len=*), intent(in)  filename,
character(len=*), intent(in)  tr_name,
real, dimension(:), intent(out), allocatable  z_edges,
integer, intent(out)  nz_out,
logical, intent(out)  has_edges,
logical, intent(inout)  use_missing,
real, intent(inout)  missing,
real, intent(in)  scale 
)
private

This subroutine reads the vertical coordinate data for a field from a NetCDF file. It also might read the missing value attribute for that same field.

Parameters
[in]filenameThe name of the file to read from.
[in]tr_nameThe name of the tracer in the file.
[out]z_edgesThe depths of the vertical edges of the tracer array
[out]nz_outThe number of vertical layers in the tracer array
[out]has_edgesIf true the values in z_edges are the edges of the tracer cells, otherwise they are the cell centers
[in,out]use_missingIf false on input, see whether the tracer has a missing value, and if so return true
[in,out]missingThe missing value, if one has been found
[in]scaleA scaling factor for z_edges into new units.

Definition at line 385 of file MOM_tracer_Z_init.F90.

385  character(len=*), intent(in) :: filename !< The name of the file to read from.
386  character(len=*), intent(in) :: tr_name !< The name of the tracer in the file.
387  real, dimension(:), allocatable, &
388  intent(out) :: z_edges !< The depths of the vertical edges of the tracer array
389  integer, intent(out) :: nz_out !< The number of vertical layers in the tracer array
390  logical, intent(out) :: has_edges !< If true the values in z_edges are the edges of the
391  !! tracer cells, otherwise they are the cell centers
392  logical, intent(inout) :: use_missing !< If false on input, see whether the tracer has a
393  !! missing value, and if so return true
394  real, intent(inout) :: missing !< The missing value, if one has been found
395  real, intent(in) :: scale !< A scaling factor for z_edges into new units.
396 
397  ! This subroutine reads the vertical coordinate data for a field from a
398  ! NetCDF file. It also might read the missing value attribute for that same field.
399  character(len=32) :: mdl
400  character(len=120) :: dim_name, edge_name, tr_msg, dim_msg
401  logical :: monotonic
402  integer :: ncid, status, intid, tr_id, layid, k
403  integer :: nz_edge, ndim, tr_dim_ids(nf90_max_var_dims)
404 
405  mdl = "MOM_tracer_Z_init read_Z_edges: "
406  tr_msg = trim(tr_name)//" in "//trim(filename)
407 
408  status = nf90_open(filename, nf90_nowrite, ncid)
409  if (status /= nf90_noerr) then
410  call mom_error(warning,mdl//" Difficulties opening "//trim(filename)//&
411  " - "//trim(nf90_strerror(status)))
412  nz_out = -1 ; return
413  endif
414 
415  status = nf90_inq_varid(ncid, tr_name, tr_id)
416  if (status /= nf90_noerr) then
417  call mom_error(warning,mdl//" Difficulties finding variable "//&
418  trim(tr_msg)//" - "//trim(nf90_strerror(status)))
419  nz_out = -1 ; status = nf90_close(ncid) ; return
420  endif
421  status = nf90_inquire_variable(ncid, tr_id, ndims=ndim, dimids=tr_dim_ids)
422  if (status /= nf90_noerr) then
423  call mom_error(warning,mdl//" cannot inquire about "//trim(tr_msg))
424  elseif ((ndim < 3) .or. (ndim > 4)) then
425  call mom_error(warning,mdl//" "//trim(tr_msg)//&
426  " has too many or too few dimensions.")
427  nz_out = -1 ; status = nf90_close(ncid) ; return
428  endif
429 
430  if (.not.use_missing) then
431  ! Try to find the missing value from the dataset.
432  status = nf90_get_att(ncid, tr_id, "missing_value", missing)
433  if (status /= nf90_noerr) use_missing = .true.
434  endif
435 
436  ! Get the axis name and length.
437  status = nf90_inquire_dimension(ncid, tr_dim_ids(3), dim_name, len=nz_out)
438  if (status /= nf90_noerr) then
439  call mom_error(warning,mdl//" cannot inquire about dimension(3) of "//&
440  trim(tr_msg))
441  endif
442 
443  dim_msg = trim(dim_name)//" in "//trim(filename)
444  status = nf90_inq_varid(ncid, dim_name, layid)
445  if (status /= nf90_noerr) then
446  call mom_error(warning,mdl//" Difficulties finding variable "//&
447  trim(dim_msg)//" - "//trim(nf90_strerror(status)))
448  nz_out = -1 ; status = nf90_close(ncid) ; return
449  endif
450  ! Find out if the Z-axis has an edges attribute
451  status = nf90_get_att(ncid, layid, "edges", edge_name)
452  if (status /= nf90_noerr) then
453  call mom_mesg(mdl//" "//trim(dim_msg)//&
454  " has no readable edges attribute - "//trim(nf90_strerror(status)))
455  has_edges = .false.
456  else
457  has_edges = .true.
458  status = nf90_inq_varid(ncid, edge_name, intid)
459  if (status /= nf90_noerr) then
460  call mom_error(warning,mdl//" Difficulties finding edge variable "//&
461  trim(edge_name)//" in "//trim(filename)//" - "//trim(nf90_strerror(status)))
462  has_edges = .false.
463  endif
464  endif
465 
466  nz_edge = nz_out ; if (has_edges) nz_edge = nz_out+1
467  allocate(z_edges(nz_edge)) ; z_edges(:) = 0.0
468 
469  if (nz_out < 1) return
470 
471  ! Read the right variable.
472  if (has_edges) then
473  dim_msg = trim(edge_name)//" in "//trim(filename)
474  status = nf90_get_var(ncid, intid, z_edges)
475  if (status /= nf90_noerr) then
476  call mom_error(warning,mdl//" Difficulties reading variable "//&
477  trim(dim_msg)//" - "//trim(nf90_strerror(status)))
478  nz_out = -1 ; status = nf90_close(ncid) ; return
479  endif
480  else
481  status = nf90_get_var(ncid, layid, z_edges)
482  if (status /= nf90_noerr) then
483  call mom_error(warning,mdl//" Difficulties reading variable "//&
484  trim(dim_msg)//" - "//trim(nf90_strerror(status)))
485  nz_out = -1 ; status = nf90_close(ncid) ; return
486  endif
487  endif
488 
489  status = nf90_close(ncid)
490  if (status /= nf90_noerr) call mom_error(warning, mdl// &
491  " Difficulties closing "//trim(filename)//" - "//trim(nf90_strerror(status)))
492 
493  ! z_edges should be montonically decreasing with our sign convention.
494  ! Change the sign sign convention if it looks like z_edges is increasing.
495  if (z_edges(1) < z_edges(2)) then
496  do k=1,nz_edge ; z_edges(k) = -z_edges(k) ; enddo
497  endif
498  ! Check that z_edges is now monotonically decreasing.
499  monotonic = .true.
500  do k=2,nz_edge ; if (z_edges(k) >= z_edges(k-1)) monotonic = .false. ; enddo
501  if (.not.monotonic) &
502  call mom_error(warning,mdl//" "//trim(dim_msg)//" is not monotonic.")
503 
504  if (scale /= 1.0) then ; do k=1,nz_edge ; z_edges(k) = scale*z_edges(k) ; enddo ; endif
505 

◆ tracer_z_init()

logical function, public mom_tracer_z_init::tracer_z_init ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(out)  tr,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
character(len=*), intent(in)  filename,
character(len=*), intent(in)  tr_name,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
real, intent(in), optional  missing_val,
real, intent(in), optional  land_val 
)

This function initializes a tracer by reading a Z-space file, returning .true. if this appears to have been successful, and false otherwise.

Returns
A return code indicating if the initialization has been successful
Parameters
[in]gThe ocean's grid structure
[in]usA dimensional unit scaling type
[out]trThe tracer to initialize
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]filenameThe name of the file to read from
[in]tr_nameThe name of the tracer in the file
[in]missing_valThe missing value for the tracer
[in]land_valA value to use to fill in land points

Definition at line 31 of file MOM_tracer_Z_init.F90.

31  logical :: tracer_z_init !< A return code indicating if the initialization has been successful
32  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
33  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
34  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
35  intent(out) :: tr !< The tracer to initialize
36  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
37  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
38  character(len=*), intent(in) :: filename !< The name of the file to read from
39  character(len=*), intent(in) :: tr_name !< The name of the tracer in the file
40 ! type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
41  real, optional, intent(in) :: missing_val !< The missing value for the tracer
42  real, optional, intent(in) :: land_val !< A value to use to fill in land points
43 
44  ! This function initializes a tracer by reading a Z-space file, returning true if this
45  ! appears to have been successful, and false otherwise.
46 !
47  integer, save :: init_calls = 0
48 ! This include declares and sets the variable "version".
49 #include "version_variable.h"
50  character(len=40) :: mdl = "MOM_tracer_Z_init" ! This module's name.
51  character(len=256) :: mesg ! Message for error messages.
52 
53  real, allocatable, dimension(:,:,:) :: &
54  tr_in ! The z-space array of tracer concentrations that is read in.
55  real, allocatable, dimension(:) :: &
56  z_edges, & ! The depths of the cell edges or cell centers (depending on
57  ! the value of has_edges) in the input z* data [Z ~> m].
58  tr_1d, & ! A copy of the input tracer concentrations in a column.
59  wt, & ! The fractional weight for each layer in the range between
60  ! k_top and k_bot, nondim.
61  z1, & ! z1 and z2 are the depths of the top and bottom limits of the part
62  z2 ! of a z-cell that contributes to a layer, relative to the cell
63  ! center and normalized by the cell thickness, nondim.
64  ! Note that -1/2 <= z1 <= z2 <= 1/2.
65  real :: e(szk_(g)+1) ! The z-star interface heights [Z ~> m].
66  real :: landval ! The tracer value to use in land points.
67  real :: sl_tr ! The normalized slope of the tracer
68  ! within the cell, in tracer units.
69  real :: htot(szi_(g)) ! The vertical sum of h [H ~> m or kg m-2].
70  real :: dilate ! The amount by which the thicknesses are dilated to
71  ! create a z-star coordinate, nondim or in m3 kg-1.
72  real :: missing ! The missing value for the tracer.
73 
74  logical :: has_edges, use_missing, zero_surface
75  character(len=80) :: loc_msg
76  integer :: k_top, k_bot, k_bot_prev, k_start
77  integer :: i, j, k, kz, is, ie, js, je, nz, nz_in
78  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
79 
80  landval = 0.0 ; if (present(land_val)) landval = land_val
81 
82  zero_surface = .false. ! Make this false for errors to be fatal.
83 
84  use_missing = .false.
85  if (present(missing_val)) then
86  use_missing = .true. ; missing = missing_val
87  endif
88 
89  ! Find out the number of input levels and read the depth of the edges,
90  ! also modifying their sign convention to be monotonically decreasing.
91  call read_z_edges(filename, tr_name, z_edges, nz_in, has_edges, use_missing, &
92  missing, scale=us%m_to_Z)
93  if (nz_in < 1) then
94  tracer_z_init = .false.
95  return
96  endif
97 
98  allocate(tr_in(g%isd:g%ied,g%jsd:g%jed,nz_in)) ; tr_in(:,:,:) = 0.0
99  allocate(tr_1d(nz_in)) ; tr_1d(:) = 0.0
100  call mom_read_data(filename, tr_name, tr_in(:,:,:), g%Domain)
101 
102  ! Fill missing values from above? Use a "close" test to avoid problems
103  ! from type-conversion rounoff.
104  if (present(missing_val)) then
105  do j=js,je ; do i=is,ie
106  if (g%mask2dT(i,j) == 0.0) then
107  tr_in(i,j,1) = landval
108  elseif (abs(tr_in(i,j,1) - missing_val) <= 1e-6*abs(missing_val)) then
109  write(loc_msg,'(f7.2," N ",f7.2," E")') g%geoLatT(i,j), g%geoLonT(i,j)
110  if (zero_surface) then
111  call mom_error(warning, "tracer_Z_init: Missing value of "// &
112  trim(tr_name)//" found in an ocean point at "//trim(loc_msg)// &
113  " in "//trim(filename) )
114  tr_in(i,j,1) = 0.0
115  else
116  call mom_error(fatal, "tracer_Z_init: Missing value of "// &
117  trim(tr_name)//" found in an ocean point at "//trim(loc_msg)// &
118  " in "//trim(filename) )
119  endif
120  endif
121  enddo ; enddo
122  do k=2,nz_in ; do j=js,je ; do i=is,ie
123  if (abs(tr_in(i,j,k) - missing_val) <= 1e-6*abs(missing_val)) &
124  tr_in(i,j,k) = tr_in(i,j,k-1)
125  enddo ; enddo ; enddo
126  endif
127 
128  allocate(wt(nz_in+1)) ; allocate(z1(nz_in+1)) ; allocate(z2(nz_in+1))
129 
130  ! This is a placeholder, and will be replaced with our full vertical
131  ! interpolation machinery when it is in place.
132  if (has_edges) then
133  do j=js,je
134  do i=is,ie ; htot(i) = 0.0 ; enddo
135  do k=1,nz ; do i=is,ie ; htot(i) = htot(i) + h(i,j,k) ; enddo ; enddo
136 
137  do i=is,ie ; if (g%mask2dT(i,j)*htot(i) > 0.0) then
138  ! Determine the z* heights of the model interfaces.
139  dilate = (g%bathyT(i,j) - 0.0) / htot(i)
140  e(nz+1) = -g%bathyT(i,j)
141  do k=nz,1,-1 ; e(k) = e(k+1) + dilate * h(i,j,k) ; enddo
142 
143  ! Create a single-column copy of tr_in. Efficiency is not an issue here.
144  do k=1,nz_in ; tr_1d(k) = tr_in(i,j,k) ; enddo
145  k_bot = 1 ; k_bot_prev = -1
146  do k=1,nz
147  if (e(k+1) > z_edges(1)) then
148  tr(i,j,k) = tr_1d(1)
149  elseif (e(k) < z_edges(nz_in+1)) then
150  tr(i,j,k) = tr_1d(nz_in)
151  else
152  k_start = k_bot ! The starting point for this search
153  call find_overlap(z_edges, e(k), e(k+1), nz_in, &
154  k_start, k_top, k_bot, wt, z1, z2)
155  kz = k_top
156  if (kz /= k_bot_prev) then
157  ! Calculate the intra-cell profile.
158  sl_tr = 0.0 ! ; cur_tr = 0.0
159  if ((kz < nz_in) .and. (kz > 1)) &
160  sl_tr = find_limited_slope(tr_1d, z_edges, kz)
161  endif
162  ! This is the piecewise linear form.
163  tr(i,j,k) = wt(kz) * (tr_1d(kz) + 0.5*sl_tr*(z2(kz) + z1(kz)))
164  ! For the piecewise parabolic form add the following...
165  ! + C1_3*cur_tr*(z2(kz)**2 + z2(kz)*z1(kz) + z1(kz)**2))
166  do kz=k_top+1,k_bot-1
167  tr(i,j,k) = tr(i,j,k) + wt(kz)*tr_1d(kz)
168  enddo
169  if (k_bot > k_top) then
170  kz = k_bot
171  ! Calculate the intra-cell profile.
172  sl_tr = 0.0 ! ; cur_tr = 0.0
173  if ((kz < nz_in) .and. (kz > 1)) &
174  sl_tr = find_limited_slope(tr_1d, z_edges, kz)
175  ! This is the piecewise linear form.
176  tr(i,j,k) = tr(i,j,k) + wt(kz) * &
177  (tr_1d(kz) + 0.5*sl_tr*(z2(kz) + z1(kz)))
178  ! For the piecewise parabolic form add the following...
179  ! + C1_3*cur_tr*(z2(kz)**2 + z2(kz)*z1(kz) + z1(kz)**2))
180  endif
181  k_bot_prev = k_bot
182 
183  ! Now handle the unlikely case where the layer partially extends
184  ! past the valid range of the input data by extrapolating using
185  ! the top or bottom value.
186  if ((e(k) > z_edges(1)) .and. (z_edges(nz_in+1) > e(k+1))) then
187  tr(i,j,k) = (((e(k) - z_edges(1)) * tr_1d(1) + &
188  (z_edges(1) - z_edges(nz_in)) * tr(i,j,k)) + &
189  (z_edges(nz_in+1) - e(k+1)) * tr_1d(nz_in)) / &
190  (e(k) - e(k+1))
191  elseif (e(k) > z_edges(1)) then
192  tr(i,j,k) = ((e(k) - z_edges(1)) * tr_1d(1) + &
193  (z_edges(1) - e(k+1)) * tr(i,j,k)) / &
194  (e(k) - e(k+1))
195  elseif (z_edges(nz_in) > e(k+1)) then
196  tr(i,j,k) = ((e(k) - z_edges(nz_in+1)) * tr(i,j,k) + &
197  (z_edges(nz_in+1) - e(k+1)) * tr_1d(nz_in)) / &
198  (e(k) - e(k+1))
199  endif
200  endif
201  enddo ! k-loop
202  else
203  do k=1,nz ; tr(i,j,k) = landval ; enddo
204  endif ; enddo ! i-loop
205  enddo ! j-loop
206  else
207  ! Without edge values, integrate a linear interpolation between cell centers.
208  do j=js,je
209  do i=is,ie ; htot(i) = 0.0 ; enddo
210  do k=1,nz ; do i=is,ie ; htot(i) = htot(i) + h(i,j,k) ; enddo ; enddo
211 
212  do i=is,ie ; if (g%mask2dT(i,j)*htot(i) > 0.0) then
213  ! Determine the z* heights of the model interfaces.
214  dilate = (g%bathyT(i,j) - 0.0) / htot(i)
215  e(nz+1) = -g%bathyT(i,j)
216  do k=nz,1,-1 ; e(k) = e(k+1) + dilate * h(i,j,k) ; enddo
217 
218  ! Create a single-column copy of tr_in. Efficiency is not an issue here.
219  do k=1,nz_in ; tr_1d(k) = tr_in(i,j,k) ; enddo
220  k_bot = 1
221  do k=1,nz
222  if (e(k+1) > z_edges(1)) then
223  tr(i,j,k) = tr_1d(1)
224  elseif (z_edges(nz_in) > e(k)) then
225  tr(i,j,k) = tr_1d(nz_in)
226  else
227  k_start = k_bot ! The starting point for this search
228  call find_overlap(z_edges, e(k), e(k+1), nz_in-1, &
229  k_start, k_top, k_bot, wt, z1, z2)
230 
231  kz = k_top
232  if (k_top < nz_in) then
233  tr(i,j,k) = wt(kz)*0.5*((tr_1d(kz) + tr_1d(kz+1)) + &
234  (tr_1d(kz+1) - tr_1d(kz))*(z2(kz)+z1(kz)))
235  else
236  tr(i,j,k) = wt(kz)*tr_1d(nz_in)
237  endif
238  do kz=k_top+1,k_bot-1
239  tr(i,j,k) = tr(i,j,k) + wt(kz)*0.5*(tr_1d(kz) + tr_1d(kz+1))
240  enddo
241  if (k_bot > k_top) then
242  kz = k_bot
243  tr(i,j,k) = tr(i,j,k) + wt(kz)*0.5*((tr_1d(kz) + tr_1d(kz+1)) + &
244  (tr_1d(kz+1) - tr_1d(kz))*(z2(kz)+z1(kz)))
245  endif
246 
247  ! Now handle the case where the layer partially extends past
248  ! the valid range of the input data.
249  if ((e(k) > z_edges(1)) .and. (z_edges(nz_in) > e(k+1))) then
250  tr(i,j,k) = (((e(k) - z_edges(1)) * tr_1d(1) + &
251  (z_edges(1) - z_edges(nz_in)) * tr(i,j,k)) + &
252  (z_edges(nz_in) - e(k+1)) * tr_1d(nz_in)) / &
253  (e(k) - e(k+1))
254  elseif (e(k) > z_edges(1)) then
255  tr(i,j,k) = ((e(k) - z_edges(1)) * tr_1d(1) + &
256  (z_edges(1) - e(k+1)) * tr(i,j,k)) / &
257  (e(k) - e(k+1))
258  elseif (z_edges(nz_in) > e(k+1)) then
259  tr(i,j,k) = ((e(k) - z_edges(nz_in)) * tr(i,j,k) + &
260  (z_edges(nz_in) - e(k+1)) * tr_1d(nz_in)) / &
261  (e(k) - e(k+1))
262  endif
263  endif
264  enddo
265  else
266  do k=1,nz ; tr(i,j,k) = landval ; enddo
267  endif ; enddo ! i-loop
268  enddo ! j-loop
269  endif
270 
271  deallocate(tr_in) ; deallocate(tr_1d) ; deallocate(z_edges)
272  deallocate(wt) ; deallocate(z1) ; deallocate(z2)
273 
274  tracer_z_init = .true.
275 

◆ tracer_z_init_array()

subroutine, public mom_tracer_z_init::tracer_z_init_array ( real, dimension(szi_(g),szj_(g),nk_data), intent(in)  tr_in,
real, dimension(nk_data+1), intent(in)  z_edges,
integer, intent(in)  nk_data,
real, dimension(szi_(g),szj_(g),nlay+1), intent(in)  e,
real, intent(in)  land_fill,
type(ocean_grid_type), intent(in)  G,
integer, intent(in)  nlay,
integer, dimension(szi_(g),szj_(g)), intent(in)  nlevs,
real, intent(in)  eps_z,
real, dimension(szi_(g),szj_(g),nlay), intent(out)  tr 
)

Layer model routine for remapping tracers from pseudo-z coordinates into layers defined by target interface positions.

Parameters
[in]gThe ocean's grid structure
[in]nk_dataThe number of levels in the input data
[in]tr_inThe z-space array of tracer concentrations that is read in.
[in]z_edgesThe depths of the cell edges in the input z* data [Z ~> m or m]
[in]nlayThe number of vertical layers in the target grid
[in]eThe depths of the target layer interfaces [Z ~> m or m]
[in]land_fillfill in data over land (1)
[in]nlevsThe number of input levels with valid data
[in]eps_zA negligibly thin layer thickness [Z ~> m].
[out]trtracers in layer space

Definition at line 282 of file MOM_tracer_Z_init.F90.

282  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
283  integer, intent(in) :: nk_data !< The number of levels in the input data
284  real, dimension(SZI_(G),SZJ_(G),nk_data), &
285  intent(in) :: tr_in !< The z-space array of tracer concentrations that is read in.
286  real, dimension(nk_data+1), intent(in) :: z_edges !< The depths of the cell edges in the input z* data
287  !! [Z ~> m or m]
288  integer, intent(in) :: nlay !< The number of vertical layers in the target grid
289  real, dimension(SZI_(G),SZJ_(G),nlay+1), &
290  intent(in) :: e !< The depths of the target layer interfaces [Z ~> m or m]
291  real, intent(in) :: land_fill !< fill in data over land (1)
292  integer, dimension(SZI_(G),SZJ_(G)), &
293  intent(in) :: nlevs !< The number of input levels with valid data
294  real, intent(in) :: eps_z !< A negligibly thin layer thickness [Z ~> m].
295  real, dimension(SZI_(G),SZJ_(G),nlay), &
296  intent(out) :: tr !< tracers in layer space
297 
298  ! Local variables
299  real, dimension(nk_data) :: tr_1d !< a copy of the input tracer concentrations in a column.
300  real, dimension(nlay+1) :: e_1d ! A 1-d column of intreface heights, in the same units as e.
301  real, dimension(nlay) :: tr_ ! A 1-d column of output tracer concentrations
302  integer :: k_top, k_bot, k_bot_prev, kstart
303  real :: sl_tr ! The tracer concentration slope times the layer thickness, in tracer units.
304  real, dimension(nk_data) :: wt !< The fractional weight for each layer in the range between z1 and z2
305  real, dimension(nk_data) :: z1, z2 ! z1 and z2 are the fractional depths of the top and bottom
306  ! limits of the part of a z-cell that contributes to a layer, relative
307  ! to the cell center and normalized by the cell thickness [nondim].
308  ! Note that -1/2 <= z1 <= z2 <= 1/2.
309  integer :: i, j, k, kz, is, ie, js, je
310 
311  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
312 
313  do j=js,je
314  i_loop: do i=is,ie
315  if (nlevs(i,j) == 0 .or. g%mask2dT(i,j) == 0.) then
316  tr(i,j,:) = land_fill
317  cycle i_loop
318  endif
319 
320  do k=1,nk_data
321  tr_1d(k) = tr_in(i,j,k)
322  enddo
323 
324  do k=1,nlay+1
325  e_1d(k) = e(i,j,k)
326  enddo
327  k_bot = 1 ; k_bot_prev = -1
328  do k=1,nlay
329  if (e_1d(k+1) > z_edges(1)) then
330  tr(i,j,k) = tr_1d(1)
331  elseif (e_1d(k) < z_edges(nlevs(i,j)+1)) then
332  tr(i,j,k) = tr_1d(nlevs(i,j))
333 
334  else
335  kstart = k_bot
336  call find_overlap(z_edges, e_1d(k), e_1d(k+1), nlevs(i,j), &
337  kstart, k_top, k_bot, wt, z1, z2)
338  kz = k_top
339  sl_tr = 0.0 ! ; cur_tr=0.0
340  if (kz /= k_bot_prev) then
341  ! Calculate the intra-cell profile.
342  if ((kz < nlevs(i,j)) .and. (kz > 1)) then
343  sl_tr = find_limited_slope(tr_1d, z_edges, kz)
344  endif
345  endif
346  if (kz > nlevs(i,j)) kz = nlevs(i,j)
347  ! This is the piecewise linear form.
348  tr(i,j,k) = wt(kz) * (tr_1d(kz) + 0.5*sl_tr*(z2(kz) + z1(kz)))
349  ! For the piecewise parabolic form add the following...
350  ! + C1_3*wt(kz) * cur_tr*(z2(kz)**2 + z2(kz)*z1(kz) + z1(kz)**2))
351  do kz=k_top+1,k_bot-1
352  tr(i,j,k) = tr(i,j,k) + wt(kz)*tr_1d(kz)
353  enddo
354 
355  if (k_bot > k_top) then
356  kz = k_bot
357  ! Calculate the intra-cell profile.
358  sl_tr = 0.0 ! ; cur_tr = 0.0
359  if ((kz < nlevs(i,j)) .and. (kz > 1)) then
360  sl_tr = find_limited_slope(tr_1d, z_edges, kz)
361  endif
362  ! This is the piecewise linear form.
363  tr(i,j,k) = tr(i,j,k) + wt(kz) * (tr_1d(kz) + 0.5*sl_tr*(z2(kz) + z1(kz)))
364  ! For the piecewise parabolic form add the following...
365  ! + C1_3*cur_tr*(z2(kz)**2 + z2(kz)*z1(kz) + z1(kz)**2))
366  endif
367  k_bot_prev = k_bot
368 
369  endif
370  enddo ! k-loop
371 
372  do k=2,nlay ! simply fill vanished layers with adjacent value
373  if (e_1d(k)-e_1d(k+1) <= eps_z) tr(i,j,k) = tr(i,j,k-1)
374  enddo
375 
376  enddo i_loop
377  enddo
378