MOM6
mom_eos::calculate_spec_vol Interface Reference

Detailed Description

Calculates specific volume of sea water from T, S and P.

Definition at line 75 of file MOM_EOS.F90.

Private functions

subroutine calc_spec_vol_scalar (T, S, pressure, specvol, EOS, spv_ref, scale)
 Calls the appropriate subroutine to calculate specific volume of sea water for scalar inputs. More...
 
subroutine calculate_spec_vol_array (T, S, pressure, specvol, start, npts, EOS, spv_ref, scale)
 Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs. More...
 
subroutine calc_spec_vol_1d (T, S, pressure, specvol, EOS, dom, spv_ref, scale)
 Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs, potentially limiting the domain of indices that are worked on. More...
 

Detailed Description

Calculates specific volume of sea water from T, S and P.

Definition at line 75 of file MOM_EOS.F90.

Functions and subroutines

◆ calc_spec_vol_1d()

subroutine mom_eos::calculate_spec_vol::calc_spec_vol_1d ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(inout)  specvol,
type(eos_type), pointer  EOS,
integer, dimension(2), intent(in), optional  dom,
real, intent(in), optional  spv_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs, potentially limiting the domain of indices that are worked on.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [R L2 T-2 ~> Pa]
[in,out]specvolIn situ specific volume [R-1 ~> m3 kg-1]
eosEquation of state structure
[in]domThe domain of indices to work on, taking into account that arrays start at 1.
[in]spv_refA reference specific volume [R-1 ~> m3 kg-1]
[in]scaleA multiplicative factor by which to scale output specific volume in combination with scaling given by US [various]

Definition at line 570 of file MOM_EOS.F90.

571  real, dimension(:), intent(in) :: T !< Potential temperature referenced to the surface [degC]
572  real, dimension(:), intent(in) :: S !< Salinity [ppt]
573  real, dimension(:), intent(in) :: pressure !< Pressure [R L2 T-2 ~> Pa]
574  real, dimension(:), intent(inout) :: specvol !< In situ specific volume [R-1 ~> m3 kg-1]
575  type(EOS_type), pointer :: EOS !< Equation of state structure
576  integer, dimension(2), optional, intent(in) :: dom !< The domain of indices to work on, taking
577  !! into account that arrays start at 1.
578  real, optional, intent(in) :: spv_ref !< A reference specific volume [R-1 ~> m3 kg-1]
579  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale
580  !! output specific volume in combination with
581  !! scaling given by US [various]
582  ! Local variables
583  real, dimension(size(specvol)) :: pres ! Pressure converted to [Pa]
584  real :: p_scale ! A factor to convert pressure to units of Pa [Pa T2 R-1 L-2 ~> 1]
585  real :: spv_unscale ! A factor to convert specific volume from R-1 to m3 kg-1 [m3 kg-1 R ~> 1]
586  real :: spv_scale ! A factor to convert specific volume from m3 kg-1 to the desired units [kg m-3 R-1 ~> 1]
587  real :: spv_reference ! spv_ref converted to [m3 kg-1]
588  integer :: i, is, ie, npts
589 
590  if (.not.associated(eos)) call mom_error(fatal, &
591  "calc_spec_vol_1d called with an unassociated EOS_type EOS.")
592 
593  if (present(dom)) then
594  is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
595  else
596  is = 1 ; ie = size(specvol) ; npts = 1 + ie - is
597  endif
598 
599  p_scale = eos%RL2_T2_to_Pa
600  spv_unscale = eos%kg_m3_to_R
601 
602  if ((p_scale == 1.0) .and. (spv_unscale == 1.0)) then
603  call calculate_spec_vol_array(t, s, pressure, specvol, is, npts, eos, spv_ref)
604  elseif (present(spv_ref)) then ! This is the same as above, but with some extra work to rescale variables.
605  do i=is,ie ; pres(i) = p_scale * pressure(i) ; enddo
606  spv_reference = spv_unscale*spv_ref
607  call calculate_spec_vol_array(t, s, pres, specvol, is, npts, eos, spv_reference)
608  else ! There is rescaling of variables, but spv_ref is not present. Passing a 0 value of spv_ref
609  ! changes answers at roundoff for some equations of state, like Wright and UNESCO.
610  do i=is,ie ; pres(i) = p_scale * pressure(i) ; enddo
611  call calculate_spec_vol_array(t, s, pres, specvol, is, npts, eos)
612  endif
613 
614  spv_scale = eos%R_to_kg_m3
615  if (present(scale)) spv_scale = spv_scale * scale
616  if (spv_scale /= 1.0) then ; do i=is,ie
617  specvol(i) = spv_scale * specvol(i)
618  enddo ; endif
619 

◆ calc_spec_vol_scalar()

subroutine mom_eos::calculate_spec_vol::calc_spec_vol_scalar ( real, intent(in)  T,
real, intent(in)  S,
real, intent(in)  pressure,
real, intent(out)  specvol,
type(eos_type), pointer  EOS,
real, intent(in), optional  spv_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate specific volume of sea water for scalar inputs.

Parameters
[in]tPotential temperature referenced to the surface [degC]
[in]sSalinity [ppt]
[in]pressurePressure [Pa] or [R L2 T-2 ~> Pa]
[out]specvolIn situ? specific volume [m3 kg-1] or [R-1 ~> m3 kg-1]
eosEquation of state structure
[in]spv_refA reference specific volume [m3 kg-1] or [R-1 m3 kg-1]
[in]scaleA multiplicative factor by which to scale specific volume in combination with scaling given by US [various]

Definition at line 532 of file MOM_EOS.F90.

533  real, intent(in) :: T !< Potential temperature referenced to the surface [degC]
534  real, intent(in) :: S !< Salinity [ppt]
535  real, intent(in) :: pressure !< Pressure [Pa] or [R L2 T-2 ~> Pa]
536  real, intent(out) :: specvol !< In situ? specific volume [m3 kg-1] or [R-1 ~> m3 kg-1]
537  type(EOS_type), pointer :: EOS !< Equation of state structure
538  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1] or [R-1 m3 kg-1]
539  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale specific
540  !! volume in combination with scaling given by US [various]
541 
542  real, dimension(1) :: Ta, Sa, pres, spv ! Rescaled single element array versions of the arguments.
543  real :: spv_reference ! spv_ref converted to [m3 kg-1]
544  real :: spv_scale ! A factor to convert specific volume from m3 kg-1 to the desired units [kg R-1 m-3 ~> 1]
545 
546  if (.not.associated(eos)) call mom_error(fatal, &
547  "calc_spec_vol_scalar called with an unassociated EOS_type EOS.")
548 
549  pres(1) = eos%RL2_T2_to_Pa*pressure
550  ta(1) = t ; sa(1) = s
551 
552  if (present(spv_ref)) then
553  spv_reference = eos%kg_m3_to_R*spv_ref
554  call calculate_spec_vol_array(ta, sa, pres, spv, 1, 1, eos, spv_reference)
555  else
556  call calculate_spec_vol_array(ta, sa, pres, spv, 1, 1, eos)
557  endif
558  specvol = spv(1)
559 
560  spv_scale = eos%R_to_kg_m3
561  if (present(scale)) spv_scale = spv_scale * scale
562  if (spv_scale /= 1.0) then
563  specvol = spv_scale * specvol
564  endif
565 

◆ calculate_spec_vol_array()

subroutine mom_eos::calculate_spec_vol::calculate_spec_vol_array ( real, dimension(:), intent(in)  T,
real, dimension(:), intent(in)  S,
real, dimension(:), intent(in)  pressure,
real, dimension(:), intent(inout)  specvol,
integer, intent(in)  start,
integer, intent(in)  npts,
type(eos_type), pointer  EOS,
real, intent(in), optional  spv_ref,
real, intent(in), optional  scale 
)
private

Calls the appropriate subroutine to calculate the specific volume of sea water for 1-D array inputs.

Parameters
[in]tpotential temperature relative to the surface [degC]
[in]ssalinity [ppt]
[in]pressurepressure [Pa]
[in,out]specvolin situ specific volume [kg m-3]
[in]startthe starting point in the arrays.
[in]nptsthe number of values to calculate.
eosEquation of state structure
[in]spv_refA reference specific volume [m3 kg-1]
[in]scaleA multiplicative factor by which to scale specific volume in combination with scaling given by US [various]

Definition at line 485 of file MOM_EOS.F90.

486  real, dimension(:), intent(in) :: T !< potential temperature relative to the surface [degC]
487  real, dimension(:), intent(in) :: S !< salinity [ppt]
488  real, dimension(:), intent(in) :: pressure !< pressure [Pa]
489  real, dimension(:), intent(inout) :: specvol !< in situ specific volume [kg m-3]
490  integer, intent(in) :: start !< the starting point in the arrays.
491  integer, intent(in) :: npts !< the number of values to calculate.
492  type(EOS_type), pointer :: EOS !< Equation of state structure
493  real, optional, intent(in) :: spv_ref !< A reference specific volume [m3 kg-1]
494  real, optional, intent(in) :: scale !< A multiplicative factor by which to scale specific
495  !! volume in combination with scaling given by US [various]
496 
497  real, dimension(size(specvol)) :: rho ! Density [kg m-3]
498  integer :: j
499 
500  if (.not.associated(eos)) call mom_error(fatal, &
501  "calculate_spec_vol_array called with an unassociated EOS_type EOS.")
502 
503  select case (eos%form_of_EOS)
504  case (eos_linear)
505  call calculate_spec_vol_linear(t, s, pressure, specvol, start, npts, &
506  eos%rho_T0_S0, eos%drho_dT, eos%drho_dS, spv_ref)
507  case (eos_unesco)
508  call calculate_spec_vol_unesco(t, s, pressure, specvol, start, npts, spv_ref)
509  case (eos_wright)
510  call calculate_spec_vol_wright(t, s, pressure, specvol, start, npts, spv_ref)
511  case (eos_teos10)
512  call calculate_spec_vol_teos10(t, s, pressure, specvol, start, npts, spv_ref)
513  case (eos_nemo)
514  call calculate_density_nemo(t, s, pressure, rho, start, npts)
515  if (present(spv_ref)) then
516  specvol(:) = 1.0 / rho(:) - spv_ref
517  else
518  specvol(:) = 1.0 / rho(:)
519  endif
520  case default
521  call mom_error(fatal, "calculate_spec_vol_array: EOS%form_of_EOS is not valid.")
522  end select
523 
524  if (present(scale)) then ; if (scale /= 1.0) then ; do j=start,start+npts-1
525  specvol(j) = scale * specvol(j)
526  enddo ; endif ; endif
527 

The documentation for this interface was generated from the following file: