Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition at line 81 of file MOM_EOS.F90.
|
| subroutine | calculate_density_derivs_scalar (T, S, pressure, drho_dT, drho_dS, EOS, scale) |
| | Calls the appropriate subroutines to calculate density derivatives by promoting a scalar to a one-element array. More...
|
| |
| subroutine | calculate_density_derivs_array (T, S, pressure, drho_dT, drho_dS, start, npts, EOS, scale) |
| | Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs. More...
|
| |
| subroutine | calculate_density_derivs_1d (T, S, pressure, drho_dT, drho_dS, EOS, dom, scale) |
| | Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs. More...
|
| |
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition at line 81 of file MOM_EOS.F90.
◆ calculate_density_derivs_1d()
| subroutine mom_eos::calculate_density_derivs::calculate_density_derivs_1d |
( |
real, dimension(:), intent(in) |
T, |
|
|
real, dimension(:), intent(in) |
S, |
|
|
real, dimension(:), intent(in) |
pressure, |
|
|
real, dimension(:), intent(inout) |
drho_dT, |
|
|
real, dimension(:), intent(inout) |
drho_dS, |
|
|
type(eos_type), pointer |
EOS, |
|
|
integer, dimension(2), intent(in), optional |
dom, |
|
|
real, intent(in), optional |
scale |
|
) |
| |
|
private |
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
- Parameters
-
| [in] | t | Potential temperature referenced to the surface [degC] |
| [in] | s | Salinity [ppt] |
| [in] | pressure | Pressure [R L2 T-2 ~> Pa] |
| [in,out] | drho_dt | The partial derivative of density with potential temperature [R degC-1 ~> kg m-3 degC-1] |
| [in,out] | drho_ds | The partial derivative of density with salinity [R degC-1 ~> kg m-3 ppt-1] |
| eos | Equation of state structure |
| [in] | dom | The domain of indices to work on, taking into account that arrays start at 1. |
| [in] | scale | A multiplicative factor by which to scale density in combination with scaling given by US [various] |
Definition at line 751 of file MOM_EOS.F90.
751 real,
dimension(:),
intent(in) :: t
752 real,
dimension(:),
intent(in) :: s
753 real,
dimension(:),
intent(in) :: pressure
754 real,
dimension(:),
intent(inout) :: drho_dt
756 real,
dimension(:),
intent(inout) :: drho_ds
758 type(eos_type),
pointer :: eos
759 integer,
dimension(2),
optional,
intent(in) :: dom
761 real,
optional,
intent(in) :: scale
764 real,
dimension(size(drho_dT)) :: pres
767 integer :: i, is, ie, npts
769 if (.not.
associated(eos))
call mom_error(fatal, &
770 "calculate_density_derivs called with an unassociated EOS_type EOS.")
772 if (
present(dom))
then 773 is = dom(1) ; ie = dom(2) ; npts = 1 + ie - is
775 is = 1 ; ie =
size(drho_dt) ; npts = 1 + ie - is
778 p_scale = eos%RL2_T2_to_Pa
780 if (p_scale == 1.0)
then 781 call calculate_density_derivs_array(t, s, pressure, drho_dt, drho_ds, is, npts, eos)
783 do i=is,ie ; pres(i) = p_scale * pressure(i) ;
enddo 784 call calculate_density_derivs_array(t, s, pres, drho_dt, drho_ds, is, npts, eos)
787 rho_scale = eos%kg_m3_to_R
788 if (
present(scale)) rho_scale = rho_scale * scale
789 if (rho_scale /= 1.0)
then ;
do i=is,ie
790 drho_dt(i) = rho_scale * drho_dt(i)
791 drho_ds(i) = rho_scale * drho_ds(i)
◆ calculate_density_derivs_array()
| subroutine mom_eos::calculate_density_derivs::calculate_density_derivs_array |
( |
real, dimension(:), intent(in) |
T, |
|
|
real, dimension(:), intent(in) |
S, |
|
|
real, dimension(:), intent(in) |
pressure, |
|
|
real, dimension(:), intent(inout) |
drho_dT, |
|
|
real, dimension(:), intent(inout) |
drho_dS, |
|
|
integer, intent(in) |
start, |
|
|
integer, intent(in) |
npts, |
|
|
type(eos_type), pointer |
EOS, |
|
|
real, intent(in), optional |
scale |
|
) |
| |
|
private |
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
- Parameters
-
| [in] | t | Potential temperature referenced to the surface [degC] |
| [in] | s | Salinity [ppt] |
| [in] | pressure | Pressure [Pa] or [R L2 T-2 ~> Pa] |
| [in,out] | drho_dt | The partial derivative of density with potential temperature [kg m-3 degC-1] or [R degC-1 ~> kg m-3 degC-1] |
| [in,out] | drho_ds | The partial derivative of density with salinity, in [kg m-3 ppt-1] or [R degC-1 ~> kg m-3 ppt-1] |
| [in] | start | Starting index within the array |
| [in] | npts | The number of values to calculate |
| eos | Equation of state structure |
| [in] | scale | A multiplicative factor by which to scale density in combination with scaling given by US [various] |
Definition at line 706 of file MOM_EOS.F90.
706 real,
dimension(:),
intent(in) :: t
707 real,
dimension(:),
intent(in) :: s
708 real,
dimension(:),
intent(in) :: pressure
709 real,
dimension(:),
intent(inout) :: drho_dt
711 real,
dimension(:),
intent(inout) :: drho_ds
713 integer,
intent(in) :: start
714 integer,
intent(in) :: npts
715 type(eos_type),
pointer :: eos
716 real,
optional,
intent(in) :: scale
722 if (.not.
associated(eos))
call mom_error(fatal, &
723 "calculate_density_derivs called with an unassociated EOS_type EOS.")
725 select case (eos%form_of_EOS)
727 call calculate_density_derivs_linear(t, s, pressure, drho_dt, drho_ds, eos%Rho_T0_S0, &
728 eos%dRho_dT, eos%dRho_dS, start, npts)
730 call calculate_density_derivs_unesco(t, s, pressure, drho_dt, drho_ds, start, npts)
732 call calculate_density_derivs_wright(t, s, pressure, drho_dt, drho_ds, start, npts)
734 call calculate_density_derivs_teos10(t, s, pressure, drho_dt, drho_ds, start, npts)
736 call calculate_density_derivs_nemo(t, s, pressure, drho_dt, drho_ds, start, npts)
738 call mom_error(fatal,
"calculate_density_derivs_array: EOS%form_of_EOS is not valid.")
741 if (
present(scale))
then ;
if (scale /= 1.0)
then ;
do j=start,start+npts-1
742 drho_dt(j) = scale * drho_dt(j)
743 drho_ds(j) = scale * drho_ds(j)
744 enddo ;
endif ;
endif
◆ calculate_density_derivs_scalar()
| subroutine mom_eos::calculate_density_derivs::calculate_density_derivs_scalar |
( |
real, intent(in) |
T, |
|
|
real, intent(in) |
S, |
|
|
real, intent(in) |
pressure, |
|
|
real, intent(out) |
drho_dT, |
|
|
real, intent(out) |
drho_dS, |
|
|
type(eos_type), pointer |
EOS, |
|
|
real, intent(in), optional |
scale |
|
) |
| |
|
private |
Calls the appropriate subroutines to calculate density derivatives by promoting a scalar to a one-element array.
- Parameters
-
| [in] | t | Potential temperature referenced to the surface [degC] |
| [in] | s | Salinity [ppt] |
| [in] | pressure | Pressure [Pa] or [R L2 T-2 ~> Pa] |
| [out] | drho_dt | The partial derivative of density with potential temperature [kg m-3 degC-1] or [R degC-1 ~> kg m-3 degC-1] |
| [out] | drho_ds | The partial derivative of density with salinity, in [kg m-3 ppt-1] or [R ppt-1 ~> kg m-3 ppt-1] |
| eos | Equation of state structure |
| [in] | scale | A multiplicative factor by which to scale density in combination with scaling given by US [various] |
Definition at line 800 of file MOM_EOS.F90.
800 real,
intent(in) :: t
801 real,
intent(in) :: s
802 real,
intent(in) :: pressure
803 real,
intent(out) :: drho_dt
805 real,
intent(out) :: drho_ds
807 type(eos_type),
pointer :: eos
808 real,
optional,
intent(in) :: scale
815 if (.not.
associated(eos))
call mom_error(fatal, &
816 "calculate_density_derivs called with an unassociated EOS_type EOS.")
818 p_scale = eos%RL2_T2_to_Pa
820 select case (eos%form_of_EOS)
822 call calculate_density_derivs_linear(t, s, p_scale*pressure, drho_dt, drho_ds, &
823 eos%Rho_T0_S0, eos%dRho_dT, eos%dRho_dS)
825 call calculate_density_derivs_wright(t, s, p_scale*pressure, drho_dt, drho_ds)
827 call calculate_density_derivs_teos10(t, s, p_scale*pressure, drho_dt, drho_ds)
829 call mom_error(fatal,
"calculate_density_derivs_scalar: EOS%form_of_EOS is not valid.")
832 rho_scale = eos%kg_m3_to_R
833 if (
present(scale)) rho_scale = rho_scale * scale
834 if (rho_scale /= 1.0)
then 835 drho_dt = rho_scale * drho_dt
836 drho_ds = rho_scale * drho_ds
The documentation for this interface was generated from the following file:
- /home/cermak/src/MOM6.devrob/src/equation_of_state/MOM_EOS.F90