///	<summary>
	///		Calculate the geopotential and temperature at the given point.
	///	</summary>
	void CalculateGeopotentialTemperature(
		const PhysicalConstants & phys,
		double dEta,
		double dXp,
		double dYp,
		double & dGeopotential,
		double & dTemperature
	) const {
		// Get some constants
		const double dG = phys.GetG();
		const double dCv = phys.GetCv();
		const double dCp = phys.GetCp();
		const double dRd = phys.GetR();
		const double dP0 = phys.GetP0();
		const double dae = phys.GetEarthRadius();
		const double df0 = 2 * phys.GetOmega() * sin(m_dRefLat);
		//const double df0 = 0.0;
		const double dbeta0 = 2 * phys.GetOmega() * cos(m_dRefLat) / dae;
		//const double dbeta0 = 0.0;
		const double dLy = m_dGDim[3] - m_dGDim[2];

		// Horizontally averaged temperature
		double dAvgTemperature =
			m_dT0 * pow(dEta, dRd * m_ddTdz / dG);

		// Horizontally averaged geopotential
		double dAvgGeopotential =
			m_dT0 * dG / m_ddTdz * 
			(1.0 - pow(dEta, dRd * m_ddTdz / dG));

		// Horizontal variation geopotential function
		double dXYGeopotential = 0.5 * m_dU0 * 
			((df0 - dbeta0 * m_dY0) * (dYp - m_dY0 - 
			m_dY0 / m_dpiC * sin(2 * m_dpiC * dYp / dLy)) + 
			0.5 * dbeta0 * 
			(dYp * dYp - dLy * dYp / m_dpiC * sin(2 * m_dpiC * dYp / dLy) - 
			0.5 * dLy * dLy / (m_dpiC * m_dpiC) * cos(2 * m_dpiC * dYp / dLy) -
			dLy * dLy / 3 - 0.5 * dLy * dLy / (m_dpiC * m_dpiC)));

		double dExpDecay = exp(-(log(dEta) / m_dbC) * (log(dEta) / m_dbC));
		double dRefProfile1 = log(dEta);
		double dRefProfile2 = 2 / (m_dbC * m_dbC) * log(dEta) * log(dEta) - 1.0;

		// Total geopotential distribution
		dGeopotential = dAvgGeopotential + dXYGeopotential*
			dRefProfile1 * dExpDecay;
		
		// Total temperature distribution
		dTemperature = dAvgTemperature + dXYGeopotential / dRd *
			dRefProfile2 * dExpDecay;
	}
	///	<summary>
	///		Evaluate the state vector at the given point.
	///	</summary>
	virtual void EvaluatePointwiseState(
		const PhysicalConstants & phys,
		const Time & time,
		double dZ,
		double dLon,
		double dLat,
		double * dState,
		double * dTracer
	) const {

		// Radius
		double dR = dZ + phys.GetEarthRadius();

		// Calculate parameters
		double dT0 = 0.5 * (ParamT0E + ParamT0P);

		double dConstA = 1.0 / ParamLapseRate;

		double dConstB = (dT0 - ParamT0P) / (dT0 * ParamT0P);

		double dConstC = 0.5 * (ParamK + 2.0)
			* (ParamT0E - ParamT0P) / (ParamT0E * ParamT0P);

		double dConstH = phys.GetR() * dT0 / phys.GetG();

		// Computed quantities
		double dScaledZ = dZ / (ParamB * dConstH);

		// Calculate tau values
		double dTau1 =
			dConstA * ParamLapseRate / dT0
				* exp(ParamLapseRate / dT0 * dZ)
			+ dConstB
				* (1.0 - 2.0 * dScaledZ * dScaledZ)
				* exp(- dScaledZ * dScaledZ);

		double dTau2 =
			dConstC * (1.0 - 2.0 * dScaledZ * dScaledZ)
				* exp(- dScaledZ * dScaledZ);

		double dIntTau1 =
			dConstA * (exp(ParamLapseRate / dT0 * dZ) - 1.0)
			+ dConstB * dZ * exp(- dScaledZ * dScaledZ);

		double dIntTau2 =
			dConstC * dZ * exp(- dScaledZ * dScaledZ);

		// Calculate utility terms
		double dRRatio;
		if (m_fDeepAtmosphere) {
			dRRatio = dR / phys.GetEarthRadius();
		} else {
			dRRatio = 1.0;
		}

		double dInteriorTerm = pow(dRRatio * cos(dLat), ParamK)
			- ParamK / (ParamK + 2.0) * pow(dRRatio * cos(dLat), ParamK + 2.0);

		// Calculate temperature
		double dTemperature = 1.0
			/ (dRRatio * dRRatio)
			/ (dTau1 - dTau2 * dInteriorTerm);

		// Calculate hydrostatic pressure
		double dPressure = phys.GetP0() * exp(
			- phys.GetG() / phys.GetR()
				* (dIntTau1 - dIntTau2 * dInteriorTerm));

		// Calculate hydrostatic density
		double dRho = dPressure / (phys.GetR() * dTemperature);

		// Velocity field
		double dInteriorTermU =
			  pow(dRRatio * cos(dLat), ParamK - 1.0)
			- pow(dRRatio * cos(dLat), ParamK + 1.0);

		double dBigU = phys.GetG() / phys.GetEarthRadius() * ParamK
			* dIntTau2 * dInteriorTermU * dTemperature;

		double dRCosLat;
		if (m_fDeepAtmosphere) {
			dRCosLat = dR * cos(dLat);
		} else {
			dRCosLat = phys.GetEarthRadius() * cos(dLat);
		}

		double dOmegaRCosLat = phys.GetOmega() * dRCosLat;

		if (dOmegaRCosLat * dOmegaRCosLat + dRCosLat * dBigU < 0.0) {
			_EXCEPTIONT("Negative discriminant detected.");
		}

		double dUlon = - dOmegaRCosLat +
			sqrt(dOmegaRCosLat * dOmegaRCosLat + dRCosLat * dBigU);

		double dUlat = 0.0;

		// Calculate velocity perturbation
		double dUlonPert;
		double dUlatPert;

		EvaluatePointwisePerturbation(
			phys, dZ, dLon, dLat,
			dUlonPert, dUlatPert);

		dUlon += dUlonPert;
		dUlat += dUlatPert;

		// Store the state
		dState[0] = dUlon;
		dState[1] = dUlat;
		dState[2] = phys.RhoThetaFromPressure(dPressure) / dRho;
		dState[3] = 0.0;
		dState[4] = dRho;

	}