Calculates the derivative of the average water density from the temperature and salinity, using a linear equation of state with the second type of averaging,
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ave_linear_eos), | intent(in) | :: | this | |||
class(scalar_field), | intent(in) | :: | temperature | A field containing the temperature of the water |
||
class(scalar_field), | intent(in) | :: | d_temperature | A field containing the derivative of the temperature of the
water, in teh same direction as |
||
class(scalar_field), | intent(in) | :: | salinity | A field containing the salinity of the water |
||
class(scalar_field), | intent(in) | :: | d_salinity | A field containing the derivative of the salinity of the
water, in the same direction as |
||
integer, | intent(in) | :: | dir | The direction in which to take the derivative |
A field containing the density of the water
function linear_water_deriv(this, temperature, d_temperature, salinity, &
d_salinity, dir) result(d_density)
!* Author: Chris MacMackin
! Date: August 2018
!
! Calculates the derivative of the average water density from the
! temperature and salinity, using a linear equation of state with
! the second type of averaging, $$ \tilde{\rho} =
! \rho_0[1-\beta_T(\tilde{\alpha}_{DT}T-T_0) +
! \beta_S(\tilde{\alpha}_{DS}S-S_0)]. $$
class(ave_linear_eos), intent(in):: this
class(scalar_field), intent(in) :: temperature
!! A field containing the temperature of the water
class(scalar_field), intent(in) :: d_temperature
!! A field containing the derivative of the temperature of the
!! water, in teh same direction as `dir`
class(scalar_field), intent(in) :: salinity
!! A field containing the salinity of the water
class(scalar_field), intent(in) :: d_salinity
!! A field containing the derivative of the salinity of the
!! water, in the same direction as `dir`
integer, intent(in) :: dir
!! The direction in which to take the derivative
class(scalar_field), pointer :: d_density
!! A field containing the density of the water
call temperature%guard_temp(); call salinity%guard_temp()
call d_temperature%guard_temp(); call d_salinity%guard_temp()
call salinity%allocate_scalar_field(d_density)
d_density = this%ref_rho*(this%a_DS_t*this%beta_s*d_salinity - &
this%a_DT_t*this%beta_t*d_temperature)
call temperature%clean_temp(); call salinity%clean_temp()
call d_temperature%clean_temp(); call d_salinity%clean_temp()
call d_density%set_temp()
end function linear_water_deriv