latentheat_all.f90

Path: physics/latentheat_all.f90
Last Update: Thu Mar 03 13:31:20 +0900 2011

    Copyright (C) GFD Dennou Club, 2004. All rights reserved.

begin

Subroutine LatentHeat

  * Developer: KITAMORI Taichi
  * Version: $Id: latentheat_all.f90,v 1.1 2011-02-23 17:26:25 yamasita Exp $
  * Tag Name: $Name:  $
  * Change History:

Overview

凝結量から凝結熱を計算する. 密度を全気相密度で評価.

Error Handling

Known Bugs

Note

Future Plans

end

Required files

Methods

Included Modules

dc_trace gridset ChemData basicset StorePotTemp

Public Instance methods

Subroutine :
xz_Masscond(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 凝結量 [kg/m^3 s]
xz_LatHeatPerMass(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 単位質量あたりの潜熱 [J/kg]
xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: エクスナー関数擾乱成分
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
: 温位擾乱成分 [K]
xz_Qcond(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
: 凝結熱による温度変化率 [K/s]

(out) 凝結熱による温度変化率

begin

Dependency

[Source]

subroutine LatentHeat(xz_Masscond, xz_LatHeatPerMass, xz_Exner, xz_PotTemp, xz_Qcond)              !(out) 凝結熱による温度変化率
                                                                 !=begin
  !==Dependency
  use dc_trace, only: BeginSub, EndSub  
  use gridset,  only: DimXMin, DimXMax, DimZMin, DimZMax
!  use physset,  only: GasR              ! 気体定数
!  use cloudset, only: SatPressB         ! Antoine の式の係数 B
  use ChemData, only: ChemData_SVapPress_AntoineB, ChemData_Cp        ! 湿潤成分の定圧比熱
  use basicset, only: xz_DensBasicZ, xz_ExnerBasicZ, xz_PotTempBasicZ, GasRDry, CpDry, CvDry, PressSfc           ! 地表面圧力
  use StorePotTemp, only: StorePotTempCond
                                                                 !=end
  !==暗黙の型宣言を禁止
  implicit none
  
  !==Input
  real(8)                :: xz_LatHeatPerMass(DimXMin:DimXMax, DimZMin:DimZMax)
                                           ! 単位質量あたりの潜熱 [J/kg]
  real(8)                :: xz_Masscond(DimXMin:DimXMax, DimZMin:DimZMax) 
                                           ! 凝結量 [kg/m^3 s]
  real(8)                :: xz_DensAll(DimXMin:DimXMax, DimZMin:DimZMax) 
                                           ! 気相密度の全量 [kg/m^3]
  real(8)                :: xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax)
                                           ! エクスナー関数擾乱成分
  real(8)                :: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) 
                                           ! 温位擾乱成分 [K]

  !==Output
  real(8), intent(out)   :: xz_Qcond(DimXMin:DimXMax, DimZMin:DimZMax) 
                                           ! 凝結熱による温度変化率 [K/s]
  
  !==Work

  call BeginSub("LatentHeat", fmt="%c", c1="Calculate potential temperature tendency bylatent heating.")

  xz_DensAll = ( PressSfc / GasRDry ) * (xz_Exner + xz_ExnerBasicZ) ** (CvDry / GasRDry) / (xz_PotTemp + xz_PotTempBasicZ)

  xz_Qcond = xz_LatHeatPerMass * xz_Masscond / (CpDry * xz_DensAll * xz_ExnerBasicZ)

  call StorePotTempCond( xz_QCond )

  call EndSub("LatentHeat")

end subroutine LatentHeat