Nodes of different colours represent the following:
Solid arrows point from a submodule to the (sub)module which it is
descended from. Dashed arrows point from a module or program unit to
modules which it uses.
Compute the potential temperature theta as a function of the temperature t, the dry air density rhod, and
the water vapor mixing ratio qv. Essentially,
.
The formulation comes from Poisson equation with equation of state plugged in and arranging
for potential temperature. This function is the exact inverse of t_of_theta_rhod_qv, which means that:
theta == theta_of_t_rhod_qv(..., t_of_theta_rhod_qv(..., theta, rhod, qv), rhod, qv).
(KCW, 2024-09-13)
Nodes of different colours represent the following:
Solid arrows point from a procedure to one which it calls. Dashed
arrows point from an interface to procedures which implement that interface.
This could include the module procedures in a generic interface or the
implementation in a submodule of an interface in a parent module.
Variables
Type
Visibility
Attributes
Name
Initial
real(kind=real64),
private
::
constant_cvd
Source Code
pure elemental function theta_of_t_rhod_qv(constant_cpd,constant_p0,constant_rd,constant_rv,t,rhod,qv)result(theta)use,intrinsic::iso_fortran_env,only:real64real(real64),intent(in)::constant_cpd,constant_p0,constant_rd,constant_rv,t,rhod,qvreal(real64)::thetareal(real64)::constant_cvd! Specific heat of dry air at constant volume.! Mayer's relation.constant_cvd=constant_cpd-constant_rd! Poisson equation with equation of state plugged in and arranging for potential temperature. For equation of state,! it can be shown that the effect of water vapor can be passed on to the temperature term entirely such that! dry air density and dry air gas constant can be used at all times. This modified "moist" temperature is! described herein:! The paragraph below equation 2.7 in doi:10.5065/1DFH-6P97.! The paragraph below equation 2 in doi:10.1175/MWR-D-11-00215.1.!! In all, solve the below equation set for $\theta$ in terms of $T$, $\rho_d$ and $q_v$:! \begin{equation*}! \begin{cases}! \theta &= T (\frac{P_0}{P})^{\frac{R_d}{C_{pd}}} \\! P &= \rho_d R_d T_m \\! T_m &= T (1 + \frac{R_v}{R_d} q_v)! \end{cases}! \end{equation*}theta=(t**(constant_cvd/constant_cpd))*&((constant_p0/(rhod*constant_rd*(1.0_real64+constant_rv/constant_rd*qv)))**&(constant_rd/constant_cpd))end function theta_of_t_rhod_qv