linear_water_density_ave2 Function

private function linear_water_density_ave2(this, temperature, salinity) result(density)

Calculates another form of the horizontally-averaged density of the water from the temperature and salinity, using a linear equation of state,

Arguments

Type IntentOptional AttributesName
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) :: salinity

A field containing the salinity of the water

Return Value class(scalar_field), pointer

A field containing the density of the water


Contents


Source Code

  function linear_water_density_ave2(this, temperature, salinity) result(density)
    !* Author: Chris MacMackin
    !  Date: August 2018
    !
    ! Calculates another form of the horizontally-averaged density of
    ! the water from the temperature and salinity, using a linear
    ! equation of state, $$ \tilde{\rho}(x) =
    ! \rho_0[1-\beta_T(\tilde{\alpha}_{DT}T(x)-T_0) +
    ! \beta_S(\tilde{\alpha}_{DS}S(x)-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)  :: salinity
      !! A field containing the salinity of the water
    class(scalar_field), pointer     :: density
      !! A field containing the density of the water
    call temperature%guard_temp(); call salinity%guard_temp()
    call salinity%allocate_scalar_field(density)
    density = this%ref_rho * (1.0_r8 - this%beta_t*(temperature - this%ref_t) &
                                     + this%beta_s*(salinity - this%ref_s))
    call temperature%clean_temp(); call salinity%clean_temp()
    call density%set_temp()
  end function linear_water_density_ave2