A linearised implementation of the equation of state which has been horizontally-integrated. The basic equation of stateis
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=r8), | private | :: | ref_rho | = | 1.0_r8 | The density for the temperature and salinity about which the equation of state was linearised, . |
|
real(kind=r8), | private | :: | ref_t | = | 0.0_r8 | The temperature about which the equation of state was linearised, . |
|
real(kind=r8), | private | :: | ref_s | = | 0.0_r8 | The salinity about which the equation of state was linearised, . |
|
real(kind=r8), | private | :: | beta_t | = | 0.0_r8 | The thermal contraction coefficient, . |
|
real(kind=r8), | private | :: | beta_s | = | 1.0_r8 | The haline contraction coefficient, . |
|
real(kind=r8), | private | :: | a_DS | = | 1.0_r8 | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction. |
|
real(kind=r8), | private | :: | a_DT | = | 1.0_r8 | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction. |
|
real(kind=r8), | private | :: | a_DS_t | = | 1.0_r8 | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction and |
|
real(kind=r8), | private | :: | a_DT_t | = | 1.0_r8 | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction and |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(in) | :: | ref_rho | The density for the temperature and salinity about which the equation of state was linearised, . |
||
real(kind=r8), | intent(in) | :: | ref_t | The temperature about which the equation of state was linearised, . |
||
real(kind=r8), | intent(in) | :: | ref_s | The salinity about which the equation of state was linearised, . |
||
real(kind=r8), | intent(in) | :: | beta_t | The thermal contraction coefficient, . |
||
real(kind=r8), | intent(in) | :: | beta_s | The haline contraction coefficient, . |
||
real(kind=r8), | intent(in), | optional | :: | a_DS | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction. Defualt value is 1. |
|
real(kind=r8), | intent(in), | optional | :: | a_DT | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction. Defualt value is 1. |
|
real(kind=r8), | intent(in), | optional | :: | a_DS_t | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction and Defualt value is 1. |
|
real(kind=r8), | intent(in), | optional | :: | a_DT_t | The shape coefficient for a horizontally-integrated model. It is defined as where and are the shapes of the variables and in the transverse direction and Defualt value is 1. |
Calculates the density of the water from the temperature and salinity, using a linear equation of state,
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) | :: | salinity | A field containing the salinity of the water |
A field containing the density of the water
Calculates one form of the horizontally-averaged density of the water from the temperature and salinity, using a linear equation of state,
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) | :: | salinity | A field containing the salinity of the water |
A field containing the density of the water
Calculates another form of the horizontally-averaged density of the water from the temperature and salinity, using a linear equation of state,
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) | :: | salinity | A field containing the salinity of the water |
A field containing the density of the water
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
Returns the haline contraction coefficient.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ave_linear_eos), | intent(in) | :: | this | |||
class(scalar_field), | intent(in) | :: | temperature | |||
class(scalar_field), | intent(in) | :: | salinity |
Returns the thermal contraction coefficient.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ave_linear_eos), | intent(in) | :: | this | |||
class(scalar_field), | intent(in) | :: | temperature | |||
class(scalar_field), | intent(in) | :: | salinity |
type, extends(equation_of_state), public :: ave_linear_eos
!* Author: Chris MacMackin
! Date: August 2018
!
! A linearised implementation of the equation of state which has
! been horizontally-integrated. The basic equation of stateis $$
! \rho(x,y) = \rho_0[1-\beta_T(T(x,y)-T_0) +
! \beta_S(S(x,y)-S_0)]. $$
!
private
real(r8) :: ref_rho = 1.0_r8
!! The density for the temperature and salinity about which the
!! equation of state was linearised, \(\rho_0\).
real(r8) :: ref_t = 0.0_r8
!! The temperature about which the equation of state was
!! linearised, \(T_0\).
real(r8) :: ref_s = 0.0_r8
!! The salinity about which the equation of state was
!! linearised, \(S_0\).
real(r8) :: beta_t = 0.0_r8
!! The thermal contraction coefficient, \(\beta_T\).
real(r8) :: beta_s = 1.0_r8
!! The haline contraction coefficient, \(\beta_S\).
real(r8) :: a_DS = 1.0_r8
!! The shape coefficient for a horizontally-integrated model. It
!! is defined as $$\alpha_{DS} = \frac{1}{y_2 - y_1}
!! \int^{y_2}_{y_1} f_{D}f_S dy, $$ where \(f_{D}(y)\) and
!! \(f_S(y)\) are the shapes of the variables \(D\) and \(S\) in
!! the transverse direction.
real(r8) :: a_DT = 1.0_r8
!! The shape coefficient for a horizontally-integrated model. It
!! is defined as $$\alpha_{DT} = \frac{1}{y_2 - y_1}
!! \int^{y_2}_{y_1} f_{D}f_T dy, $$ where \(f_{D}(y)\) and
!! \(f_T(y)\) are the shapes of the variables \(D\) and \(T\) in
!! the transverse direction.
real(r8) :: a_DS_t = 1.0_r8
!! The shape coefficient for a horizontally-integrated model. It
!! is defined as $$\tilde{\alpha}_{DS} =
!! \frac{1}{\alpha_{D^2}(y_2 - y_1)} \int^{y_2}_{y_1} f^2_{D}f_S
!! dy, $$ where \(f_{D}(y)\) and \(f_S(y)\) are the shapes of
!! the variables \(D\) and \(S\) in the transverse direction and
!! $$\alpha_{D^2} = \frac{1}{y_2 -
!! y_1}\int^{y_2}_{y_1}f_D^2dy.$$
real(r8) :: a_DT_t = 1.0_r8
!! The shape coefficient for a horizontally-integrated model. It
!! is defined as $$\tilde{\alpha}_{DT} =
!! \frac{1}{\alpha_{D^2}(y_2 - y_1)} \int^{y_2}_{y_1} f^2_{D}f_T
!! dy, $$ where \(f_{D}(y)\) and \(f_T(y)\) are the shapes of
!! the variables \(D\) and \(T\) in the transverse direction and
!! $$\alpha_{D^2} = \frac{1}{y_2 -
!! y_1}\int^{y_2}_{y_1}f_D^2dy.$$
contains
procedure :: water_density => linear_water_density
procedure :: water_density_ave1 => linear_water_density_ave1
procedure :: water_density_ave2 => linear_water_density_ave2
procedure :: water_density_derivative => linear_water_deriv
procedure :: haline_contraction => linear_haline_contraction
procedure :: thermal_contraction => linear_thermal_contraction
end type ave_linear_eos