| Path: | env/basicenv_smooth2.f90 |
| Last Update: | Mon Sep 10 10:53:27 +0900 2012 |
| Authors: | SUGIYAMA Koichiro, ODAKA Masatsugu |
| Version: | $Id: basicenv_smooth2.f90,v 1.1 2012-09-10 01:53:27 yamasita Exp $ |
| Tag Name: | $Name: arare4-20120911 $ |
| Copyright: | Copyright (C) GFD Dennou Club, 2006. All rights reserved. |
| License: | See COPYRIGHT |
デフォルトの基本場を設定するための変数参照型モジュール
* BasicEnvFile_init: 基本場の値を netCDF ファイルから取得 * BasicEnvCalc_Init: 基本場の情報を Namelist から取得して値を計算
| Subroutine : |
デフォルトの基本場を設定するためのサブルーチン. 基本場を計算し, BasicSet モジュールの値を初期化する.
コンパイルの順序の問題から, 基本場の値(hogeBasicZ な変数)を 計算する部分をBasicSet モジュールから切り離している. ECCM 始め, BasicSet 自体に依存するが hogeBasicZ は use しない 外部サブルーチンを利用するためである.
subroutine BasicEnv()
!
!デフォルトの基本場を設定するためのサブルーチン.
!基本場を計算し, BasicSet モジュールの値を初期化する.
!
!コンパイルの順序の問題から, 基本場の値(hogeBasicZ な変数)を
!計算する部分をBasicSet モジュールから切り離している.
!ECCM 始め, BasicSet 自体に依存するが hogeBasicZ は use しない
!外部サブルーチンを利用するためである.
!
!モジュール読み込み
use dc_message, only: MessageNotify
use gridset, only: DimXMin, RegXMin, DimXMax, RegXMax, DimZMin, RegZMin, DimZMax, RegZMax, SpcNum, s_Z, DelZ !Z 方向の格子点間隔
use basicset, only: BasicSetArray_Init, PressBasis, GasRDry, CpDry, CvDry, MolWtDry, Grav, SpcWetMolFr, MolWtWet, EnvType, Tropopause, GasRUniv, Humidity, TempStrat, Dhight, TempSfc, PressSfc
use Boundary, only: BoundaryXCyc_xz, BoundaryZSym_xz, BoundaryXCyc_xza, BoundaryZSym_xza !
use ECCM, only: ECCM_Dry, ECCM_Wet
use cloudset, only: SatRtWetAdia
use ChemData, only: ChemData_SVapPress_AntoineA, ChemData_SVapPress_AntoineB
!暗黙の型宣言禁止
implicit none
!変数の定義
real(8) :: xz_DensBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_PressBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_ExnerBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_TempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_PotTempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_VelSoundBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xza_MixRtBasicZ(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
real(8) :: xz_EffMolWtBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: z_TempBasicZ(DimZMin:DimZMax)
real(8) :: z_PressBasicZ(DimZMin:DimZMax)
real(8) :: MolFrIni(SpcNum)
real(8) :: xza_MolFr(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
real(8) :: za_MolFr(DimZMin:DimZMax, SpcNum)
real(8) :: xza_MixRtDivMolWt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
real(8) :: z_DTempDZ(DimZMin:DimZMax)
real(8) :: z_MolWtMean(DimZMin:DimZMax)
real(8) :: DTempDZ(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: Weight1(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: Weight2(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: LapseRateRatio ! 乾燥断熱減率に対する下層での温度減率の比
! 乾燥断熱よりも小さく設定する.
real(8) :: xz_TempAdia(DimXMin:DimXMax, DimZMin:DimZMax)
! 乾燥断熱線に沿った温度
real(8) :: xz_TempSat(DimXMin:DimXMax, DimZMin:DimZMax)
! 凝結線に沿った温度
real(8) :: xz_TempIso(DimXMin:DimXMax, DimZMin:DimZMax)
! 等温線に沿った温度
real(8) :: xz_PressAdia(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_PressSat(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_PressIso(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_TempWork(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_PressWork(DimXMin:DimXMax, DimZMin:DimZMax)
real(8) :: xz_Z(DimXMin:DimXMax, DimZMin:DimZMax)
! 2D 座標
real(8) :: RKWork1(DimXMin:DimXMax, DimZMin:DimZMax)
! ルンゲ・クッタ法で計算を行なう際に用いる作業変数
real(8) :: RKWork2(DimXMin:DimXMax, DimZMin:DimZMax)
! ルンゲ・クッタ法で計算を行なう際に用いる作業変数
real(8) :: RKWork3(DimXMin:DimXMax, DimZMin:DimZMax)
! ルンゲ・クッタ法で計算を行なう際に用いる作業変数
real(8) :: RKWork4(DimXMin:DimXMax, DimZMin:DimZMax)
! ルンゲ・クッタ法で計算を行なう際に用いる作業変数
real(8) :: Temp_0, Temp_1, Press_0, Press_1
! 乾燥断熱線と湿潤断熱線とが交わる高度を反復法で求める
! 際に用いる作業変数
real(8) :: Work ! 作業変数
real(8) :: LCL ! 下層の温度分布と湿潤断熱線が交わる高度
real(8) :: LTP ! 湿潤断熱線と等温線が交わる高度
real(8) :: Gamma_1, Gamma_2, Gamma_3
! 3 つの一次関数で温度分布を接続する際の各温度減率
real(8) :: HTrans
! 接続を開始する高度
real(8) :: LCLFunc
! LCL を求める為の作業関数
real(8) :: LCLFuncDeriv
! LCL を求める為の作業関数の導関数
real(8) :: Height_1
real(8) :: Height_2
! LCL を求める為の作業変数
! real(8) :: xz_LCL(DimXMin:DimXMax, DimZMin:DimZMax)
! ! 下層の温度分布と湿潤断熱線が交わる高度
! real(8) :: xz_LTP(DimXMin:DimXMax, DimZMin:DimZMax)
! ! 下層の温度分布と湿潤断熱線が交わる高度
integer :: i, k, s
!---------------------------------------------------------------
! 配列の初期化
!---------------------------------------------------------------
xz_PressBasicZ = 0.0d0
xz_ExnerBasicZ = 0.0d0
xz_TempBasicZ = 0.0d0
xz_PotTempBasicZ = 0.0d0
xz_VelSoundBasicZ = 0.0d0
xza_MixRtBasicZ = 0.0d0
xz_EffMolWtBasicZ = 0.0d0
! z_TempBasicZ = 0.0d0
z_TempBasicZ = TempSfc
z_PressBasicZ = 0.0d0
za_MolFr = 0.0d0
! 座標の初期化
do k = DimZMin, DimZMax
xz_Z(:,k) = s_Z(k)
end do
!---------------------------------------------------------------
! EnvType を元に, 温度, 圧力, 組成を決める
!---------------------------------------------------------------
MolFrIni = SpcWetMolFr(1:SpcNum)
!場合分け. Dry なら乾燥断熱的に, Wet なら湿潤断熱的な初期場を決める
select case(EnvType)
case("Dry")
call ECCM_Dry( MolFrIni, Humidity, z_TempBasicZ, z_PressBasicZ, z_MolWtMean, za_MolFr )
case("Wet")
call ECCM_Wet( MolFrIni, Humidity, z_TempBasicZ, z_PressBasicZ, z_MolWtMean, za_MolFr )
end select
! 北守修論計算の基本場設定(極冠周縁での典型的温度プロファイル)を改良したもの.
! * 温度勾配の不連続を解消する為, 重み付き関数でなめらかに接続.
! * 基本場の温度分布が超断熱, 過飽和とならないよう調整.
if (TempSfc <= 300.0d0 ) then
! 乾燥断熱減率に対する下層の温度減率の比
LapseRateRatio = 1.0d0
! xz_TempAdia = TempSfc - Grav * xz_Z / CpDry
xz_TempAdia = TempSfc - LapseRateRatio * Grav * xz_Z / CpDry
xz_TempIso = TempStrat
!--- 乾燥断熱線, 湿潤断熱線, 等温線が交わる高度を計算し,
!--- 各領域で成り立つ式を用いて温度, 圧力を計算
!--- 乾燥断熱線と湿潤断熱線が交わる高度(LCL)を反復法で計算
Press_0 = PressSfc
Temp_0 = TempSfc
do
Temp_1 = ChemData_SVapPress_AntoineB(12) / (ChemData_SVapPress_AntoineA(12) - dlog(Press_0/SatRtWetAdia))
Press_1 = PressSfc * (Temp_1/TempSfc)**(CpDry / (LapseRateRatio * GasRDry))
if (abs(Temp_1 - Temp_0) < epsilon(0.0d0)) then
LCL = (TempSfc * CpDry) / (Grav * LapseRateRatio) * (1.0d0 - (Press_1/PressSfc)**( GasRDry/CpDry ))
exit
else
Temp_0 = Temp_1
Press_0 = Press_1
end if
end do
!--- 湿潤断熱線と等温線が交わる高度(LTP)を計算
LTP = LCL + GasRDry * ChemData_SVapPress_AntoineB(12) / Grav * dlog(Temp_1/xz_TempIso(1,1))
!--- LCL, LTP の値を表示(後での Newton 法の初期値に用いる)
write(*,*) 'LCL', LCL
write(*,*) 'LTP', LTP
!--- 多段階的に接続するのに用いる 1 次関数の温度減率を計算
Gamma_1 = 0.75d0 * Grav / CpDry + 0.25d0 * Grav / GasRDry / ChemData_SVapPress_AntoineB(12) * Temp_1
Gamma_2 = 0.50d0 * Grav / CpDry + 0.50d0 * Grav / GasRDry / ChemData_SVapPress_AntoineB(12) * Temp_1
Gamma_3 = 0.25d0 * Grav / CpDry + 0.75d0 * Grav / GasRDry / ChemData_SVapPress_AntoineB(12) * Temp_1
write(*,*) 'Gamma1', Gamma_1
write(*,*) 'Gamma2', Gamma_2
write(*,*) 'Gamma3', Gamma_3
!--- 接続を開始する高度
HTrans = LCL - 3.0 * DelZ
!--- なめらかに接続したときの LCL を Newton 法で求める
Height_1 = LCL
do
LCLFunc = dlog(PressSfc) + CpDry / GasRDry * dlog(1.0d0 - Grav * HTrans / CpDry / TempSfc) + Grav / GasRDry / Gamma_1 * ( dlog( - Height_1 / 3.0d0 + TempSfc / Gamma_1 + (1.0d0/3.0d0 - Grav / CpDry / Gamma_1) * HTrans ) - dlog( TempSfc / Gamma_1 - Grav / CpDry / Gamma_1 * HTrans) ) + Grav / GasRDry / Gamma_2 * ( dlog( - (1.0d0 + Gamma_1 / Gamma_2) * Height_1 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 + 1.0d0/3.0d0 ) * HTrans ) - dlog( - Gamma_1 * Height_1 / Gamma_2 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 ) * HTrans ) ) + Grav / GasRDry / Gamma_3 * ( dlog( - (1.0d0 + Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) * Height_1 / 3.0d0 + TempSfc / Gamma_3 + ( - Grav / CpDry / Gamma_3 + Gamma_1 / Gamma_3 / 3.0d0 + Gamma_2 / Gamma_3 / 3.0d0 + 1.0d0/3.0d0 ) * HTrans ) - dlog( - (Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) * Height_1 / 3.0d0 + TempSfc / Gamma_3 + ( - Grav / CpDry / Gamma_3 + Gamma_1 / Gamma_3 / 3.0d0 + Gamma_2 / Gamma_3 / 3.0d0) * HTrans ) ) - ChemData_SVapPress_AntoineA(12) + ChemData_SVapPress_AntoineB(12) / ( - ( Gamma_1 + Gamma_2 + Gamma_3 ) * Height_1 / 3.0d0 + TempSfc + ( - Grav / CpDry + ( Gamma_1 + Gamma_2 + Gamma_3 ) / 3.0d0) * HTrans )
write(*,*) 'LCLFunc', LCLFunc
LCLFuncDeriv = - Grav / GasRDry / Gamma_1 / 3.0d0 / ( - Height_1 / 3.0d0 + TempSfc / Gamma_1 + (1.0d0/3.0d0 - Grav / CpDry / Gamma_1) * HTrans ) + Grav / GasRDry / Gamma_2 * ( - (1.0d0 + Gamma_1 / Gamma_2) / 3.0d0 / ( - (1.0d0 + Gamma_1 / Gamma_2) * Height_1 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 + 1.0d0/3.0d0 ) * HTrans ) + Gamma_1 / Gamma_2 / 3.0d0 / ( - Gamma_1 * Height_1 / Gamma_2 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 ) * HTrans ) ) + Grav / GasRDry / Gamma_3 * ( - (1.0d0 + Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) / 3.0d0 / ( - (1.0d0 + Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) * Height_1 / 3.0d0 + TempSfc / Gamma_3 + ( - Grav / CpDry / Gamma_3 + Gamma_1 / Gamma_3 / 3.0d0 + Gamma_2 / Gamma_3 / 3.0d0 + 1.0d0/3.0d0 * HTrans ) ) + (Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) / 3.0d0 / ( - (Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) * Height_1 / 3.0d0 + TempSfc / Gamma_3 + ( - Grav / CpDry / Gamma_3 + Gamma_1 / Gamma_3 / 3.0d0 + Gamma_2 / Gamma_3 / 3.0d0) * HTrans ) ) + ( Gamma_1 + Gamma_2 + Gamma_3 ) * ChemData_SVapPress_AntoineB(12) / 3.0d0 / ( - ( Gamma_1 + Gamma_2 + Gamma_3 ) * Height_1 / 3.0d0 + TempSfc + ( - Grav / CpDry + ( Gamma_1 + Gamma_2 + Gamma_3 ) / 3.0d0) * HTrans)**2
Height_2 = Height_1 - LCLFunc / LCLFuncDeriv
if (abs(Height_2 - Height_1) < epsilon(0.0e0)) then
LCL = Height_2
Temp_1 = - ( Gamma_1 + Gamma_2 + Gamma_3 ) * LCL / 3.0d0 + TempSfc + ( - Grav / CpDry + ( Gamma_1 + Gamma_2 + Gamma_3 ) / 3.0d0) * HTrans
write(*,*) 'LCL', LCL
write(*,*) 'TLCL', Temp_1
exit
else
Height_1 = Height_2
write(*,*) 'Height1', Height_1
end if
end do
!--- 湿潤断熱線と等温線が交わる高度(LTP)を再計算
LTP = LCL + GasRDry * ChemData_SVapPress_AntoineB(12) / Grav * dlog(Temp_1/xz_TempIso(1,1))
write(*,*) 'LTP', LTP
do k = DimZMin, DimZMax
if (s_Z(k) <= HTrans ) then
xz_TempBasicZ(:,k) = xz_TempAdia(:,k)
else if ( (s_Z(k) >= HTrans) .and. (s_Z(k) < LCL/3.0d0 + 2.0d0*HTrans/3.0d0)) then
xz_TempBasicZ(:,k) = - Gamma_1 * s_Z(k) + TempSfc + (Gamma_1 - Grav / CpDry) * HTrans
else if ( (s_Z(k) >= LCL/3.0d0 + 2.0d0*HTrans/3.0d0) .and. (s_Z(k) < 2.0d0*LCL/3.0d0 + HTrans/3.0d0)) then
xz_TempBasicZ(:,k) = - Gamma_2 * s_Z(k) + (Gamma_2 - Gamma_1) * LCL / 3.0d0 + TempSfc + (- Grav / CpDry + Gamma_1 / 3.0d0 + 2.0d0 * Gamma_2 / 3.0d0) * HTrans
else if ( (s_Z(k) >= 2.0d0*LCL/3.0d0 + HTrans/3.0d0) .and. (s_Z(k) < LCL)) then
xz_TempBasicZ(:,k) = - Gamma_3 * s_Z(k) + (2.0d0 * Gamma_3 - Gamma_2 - Gamma_1 ) * LCL / 3.0d0 + TempSfc + (- Grav / CpDry + Gamma_1 / 3.0d0 + Gamma_2 / 3.0d0 + Gamma_3 / 3.0d0) * HTrans
else if ((s_Z(k) >= LCL) .and. (s_Z(k) < LTP)) then
xz_TempBasicZ(:,k) = Temp_1 * dexp(-Grav * (s_Z(k) - LCL) /(GasRDry*ChemData_SVapPress_AntoineB(12)))
else
xz_TempBasicZ(:,k) = xz_TempIso(:,k)
end if
end do
else
stop 'Tha value of TempSfc is invalid.'
end if
! 圧力分布は解析的に求めたものを与える.
do k = DimZMin, DimZMax
if (s_Z(k) <= HTrans ) then
xz_PressBasicZ(:,k) = PressSfc * (1.0d0 - Grav * s_Z(k) / CpDry / TempSfc)**(CpDry/GasRDry)
else if ( (s_Z(k) >= HTrans) .and. (s_Z(k) <= LCL/3.0d0 + 2.0d0*HTrans/3.0d0)) then
xz_PressBasicZ(:,k) = PressSfc * (1.0d0 - Grav * HTrans / CpDry / TempSfc)**(CpDry/GasRDry) * ( (- s_Z(k) + (TempSfc + (Gamma_1 - Grav/CpDry) * HTrans)/Gamma_1 ) / ( TempSfc / Gamma_1 - Grav / CpDry / Gamma_1 * HTrans ) )**(Grav / GasRDry / Gamma_1)
else if ( (s_Z(k) >= LCL/3.0d0 + 2.0d0*HTrans/3.0d0) .and. (s_Z(k) <= 2.0d0*LCL/3.0d0 + HTrans/3.0d0)) then
xz_PressBasicZ(:,k) = PressSfc * (1.0d0 - Grav * HTrans / CpDry / TempSfc)**(CpDry/GasRDry) * ( (- LCL / 3.0d0 + TempSfc / Gamma_1 + (1.0d0/3.0d0 - Grav / CpDry / Gamma_1) * HTrans) / ( TempSfc / Gamma_1 - Grav / CpDry / Gamma_1 * HTrans ) )**(Grav / GasRDry / Gamma_1) * ( (- s_Z(k) + ( (Gamma_2 - Gamma_1)*LCL/3.0d0 + TempSfc + ( - Grav/CpDry + Gamma_1/3.0d0 + 2.0d0 * Gamma_2/3.0d0 )*HTrans ) /Gamma_2) / ( - Gamma_1 * LCL / Gamma_2 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 ) * HTrans ) )**(Grav / GasRDry / Gamma_2)
else if ( (s_Z(k) >= 2.0d0*LCL/3.0d0 + HTrans/3.0d0) .and. (s_Z(k) <= LCL)) then
xz_PressBasicZ(:,k) = PressSfc * (1.0d0 - Grav * HTrans / CpDry / TempSfc)**(CpDry/GasRDry) * ( (- LCL / 3.0d0 + TempSfc / Gamma_1 + (1.0d0/3.0d0 - Grav / CpDry / Gamma_1) * HTrans) / (TempSfc / Gamma_1 - Grav / CpDry / Gamma_1 * HTrans) )**(Grav / GasRDry / Gamma_1) * ( ( - (1.0d0 + Gamma_1 / Gamma_2) * LCL / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 + 1.0d0/3.0d0 ) * HTrans) / ( - Gamma_1 * LCL / Gamma_2 / 3.0d0 + TempSfc / Gamma_2 + ( Gamma_1 / Gamma_2 / 3.0d0 - Grav / CpDry / Gamma_2 ) * HTrans ) )**(Grav / GasRDry / Gamma_2) * ( (- s_Z(k) + ( (2.0d0*Gamma_3 - Gamma_2 - Gamma_1)*LCL/3.0d0 + TempSfc + ( - Grav/CpDry + (Gamma_1 + Gamma_2 + Gamma_3)/3.0d0 )*HTrans)/Gamma_3) / ( - (Gamma_2 / Gamma_3 + Gamma_1 / Gamma_3) * LCL / 3.0d0 + TempSfc / Gamma_3 + ( - Grav / CpDry / Gamma_3 + Gamma_1 / Gamma_3 / 3.0d0 + Gamma_2 / Gamma_3 / 3.0d0) * HTrans) )**(Grav / GasRDry / Gamma_3)
else if ((s_Z(k) >= LCL) .and. (s_Z(k) <= LTP)) then
xz_PressBasicZ(:,k) = dexp( ChemData_SVapPress_AntoineA(12) - ChemData_SVapPress_AntoineB(12) / xz_TempBasicZ(:,k))
else
xz_PressBasicZ(:,k) = dexp( ChemData_SVapPress_AntoineA(12) - ChemData_SVapPress_AntoineB(12) / xz_TempIso(1,1)) * dexp( - Grav / GasRDry / xz_TempIso(1,1) * (s_Z(k) - LTP))
end if
end do
!確認のため出力
call MessageNotify( "M", "BasicEnv", "Basic State Atmospheric Profiles." )
do k = RegZMin+1, DimZMax-1
write(*,*) "temp", k, s_Z(k), xz_TempBasicZ(1,k), xz_PressBasicZ(1,k)
end do
! 2 次元配列に格納
! do i = DimXMin, DimXMax
! xz_TempBasicZ(i,:) = z_TempBasicZ
! xz_PressBasicZ(i,:) = z_PressBasicZ
! end do
!境界条件
call BoundaryXCyc_xz( xz_TempBasicZ )
call BoundaryZSym_xz( xz_TempBasicZ )
call BoundaryXCyc_xz( xz_PressBasicZ )
call BoundaryZSym_xz( xz_PressBasicZ )
!---------------------------------------------------------------
! 混合比
!---------------------------------------------------------------
!水平方向には一様
do i = DimXMin, DimXMax
xza_MolFr(i,:,:) = za_MolFr
end do
!気相のモル比を混合比に変換
do s = 1, SpcNum
xza_MixRtBasicZ(:,:,s) = xza_MolFr(:,:,s) * MolWtWet(s) / MolWtDry
end do
! !値が小さくなりすぎないように最低値を与える
! where (xza_MixRtBasicZ <= 1.0d-20 )
! xza_MixRtBasicZ = 1.0d-20
! end where
!境界条件
call BoundaryXCyc_xza( xza_MixRtBasicZ )
call BoundaryZSym_xza( xza_MixRtBasicZ )
!---------------------------------------------------------------
! 分子量の効果
!---------------------------------------------------------------
do s = 1, SpcNum
xza_MixRtDivMolWt(:,:,s) = xza_MixRtBasicZ(:,:,s) / MolWtWet(s)
end do
xz_EffMolWtBasicZ = (1.0d0 + sum(xza_MixRtBasicZ,3) ) / ( MolWtDry * ((1.0d0 / MolWtDry) + sum(xza_MixRtDivMolWt,3)) )
!境界条件
call BoundaryXCyc_xz( xz_EffMolWtBasicZ )
call BoundaryZSym_xz( xz_EffMolWtBasicZ )
!---------------------------------------------------------------
! 温位
!---------------------------------------------------------------
xz_PotTempBasicZ = xz_TempBasicZ * (PressBasis / xz_PressBasicZ) ** (GasRDry / CpDry)
!境界条件
call BoundaryXCyc_xz( xz_PotTempBasicZ )
call BoundaryZSym_xz( xz_PotTempBasicZ )
!---------------------------------------------------------------
! エクスナー関数
!---------------------------------------------------------------
xz_ExnerBasicZ = xz_TempBasicZ / xz_PotTempBasicZ
!境界条件
call BoundaryXCyc_xz( xz_ExnerBasicZ )
call BoundaryZSym_xz( xz_ExnerBasicZ )
!---------------------------------------------------------------
! 密度
!---------------------------------------------------------------
xz_DensBasicZ = PressBasis * (xz_ExnerBasicZ ** (CvDry / GasRDry)) / (GasRDry * xz_PotTempBasicZ / xz_EffMolWtBasicZ)
!境界条件
call BoundaryXCyc_xz( xz_DensBasicZ )
call BoundaryZSym_xz( xz_DensBasicZ )
!---------------------------------------------------------------
! 音速
!---------------------------------------------------------------
xz_VelSoundBasicZ = sqrt ( CpDry * GasRDry * xz_ExnerBasicZ * xz_PotTempBasicZ / (CvDry * xz_EffMolWtBasicZ) )
!境界条件
call BoundaryXCyc_xz( xz_VelSoundBasicZ )
call BoundaryZSym_xz( xz_VelSoundBasicZ )
!----------------------------------------------------------
! BasicSet モジュールに値を設定
!----------------------------------------------------------
call BasicSetArray_Init( xz_PressBasicZ, xz_ExnerBasicZ, xz_TempBasicZ, xz_PotTempBasicZ, xz_DensBasicZ, xz_VelSoundBasicZ, xza_MixRtBasicZ, xz_EffMolWtBasicZ )
end subroutine BasicEnv