| Class | vdiffusion_jma_my_wrapper |
| In: |
vdiffusion/vdiffusion_jma_my_wrapper.F90
|
Note that Japanese and English are described in parallel.
鉛直拡散フラックスを計算します.
Vertical diffusion flux is calculated.
| VDiffusion : | 鉛直拡散フラックスの計算 |
| VDiffusionOutPut : | フラックスの出力 |
| ———— : | ———— |
| VDiffusion : | Calculate vertical diffusion fluxes |
| VDiffusionOutPut : | Output fluxes |
| Subroutine : |
vdiffusion_my モジュールの初期化を行います. NAMELIST#vdiffusion_my_nml の読み込みはこの手続きで行われます.
"vdiffusion_my" module is initialized. "NAMELIST#vdiffusion_my_nml" is loaded in this procedure.
subroutine VDiffusionJMAInit
!
! vdiffusion_my モジュールの初期化を行います.
! NAMELIST#vdiffusion_my_nml の読み込みはこの手続きで行われます.
!
! "vdiffusion_my" module is initialized.
! "NAMELIST#vdiffusion_my_nml" is loaded in this procedure.
!
! モジュール引用 ; USE statements
!
! NAMELIST ファイル入力に関するユーティリティ
! Utilities for NAMELIST file input
!
!!$ use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid
! ファイル入出力補助
! File I/O support
!
use dc_iounit, only: FileOpen
! 種別型パラメタ
! Kind type parameter
!
use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output
! 文字列操作
! Character handling
!
use dc_string, only: StoA
#ifdef JMAVDIFF
! 時刻管理
! Time control
!
use timeset, only: DelTime ! $ \Delta t $ [s]
#endif
! ヒストリデータ出力
! History data output
!
use gtool_historyauto, only: HistoryAutoAddVariable
! 補助的な変数を計算するサブルーチン・関数群
! Subroutines and functions for calculating auxiliary variables
!
use auxiliary, only : RefPressForPotTemp => RefPress, AuxVarsInit
#ifdef JMAVDIFF
!
! JMA library
!
use pp_phys_const, only: pp_phys_const_set
use pbl_mym_option_symbol, only: mymodel25, mymodel3
use pbl_grid, only: pbl_grid_ini
use pbl_const, only: pbl_const_ini
use pbl_parm, only: pbl_parm_ini
use pbl_mym_option, only: pbl_mym_option_ini
use pbl_mym_parm, only: pbl_mym_parm_ini
use pbl_mym_const, only: pbl_mym_const_ini
#endif
! 座標データ設定
! Axes data settings
!
use axesset, only: AxnameX, AxnameY, AxnameZ, AxnameR, AxnameT
! 宣言文 ; Declaration statements
!
implicit none
#ifdef JMAVDIFF
real(DP) :: PressRef
real(DP) :: TimeStep
real(DP) :: SurfEmis
#endif
!!$ integer:: unit_nml ! NAMELIST ファイルオープン用装置番号.
!!$ ! Unit number for NAMELIST file open
!!$ integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT.
!!$ ! IOSTAT of NAMELIST read
! NAMELIST 変数群
! NAMELIST group name
!
!!$ namelist /vdiffusion_my_nml/ &
!!$ & FlagConstDiffCoef, &
!!$ & ConstDiffCoefM, ConstDiffCoefH, &
!!$!
!!$ & SquareVelMin, BulkRiNumMin, &
!!$!
!!$ & MixLengthMax, ShMin, SmMin, &
!!$ & VelDiffCoefMin, TempDiffCoefMin, &
!!$ & VelDiffCoefMax, TempDiffCoefMax, &
!!$!
!!$ & MYConstA1, MYConstB1, MYConstA2, MYConstB2, MYConstC1
!
! デフォルト値については初期化手続 "vdiffusion_my#VDiffInit"
! のソースコードを参照のこと.
!
! Refer to source codes in the initialization procedure
! "vdiffusion_my#VDiffInit" for the default values.
!
! 実行文 ; Executable statement
!
if ( vdiffusion_jma_my_wrap_inited ) return
! デフォルト値の設定
! Default values settings
!
!!$ FlagConstDiffCoef = .false.
!!$ ConstDiffCoefM = 0.0_DP
!!$ ConstDiffCoefH = 0.0_DP
!!$
!!$ SquareVelMin = 0.1_DP
!!$ BulkRiNumMin = - 100.0_DP
!!$
!!$ MixLengthMax = 300.0_DP
!!$ ShMin = 0.0_DP
!!$ SmMin = 0.0_DP
!!$ VelDiffCoefMin = 0.1_DP
!!$ TempDiffCoefMin = 0.1_DP
!!$ VelDiffCoefMax = 10000.0_DP
!!$ TempDiffCoefMax = 10000.0_DP
!!$
!!$ ! Parameters proposed by Mellor and Yamada (1982).
!!$ !
!!$ MYConstA1 = 0.92_DP
!!$ MYConstB1 = 16.6_DP
!!$ MYConstA2 = 0.74_DP
!!$ MYConstB2 = 10.1_DP
!!$ MYConstC1 = 0.08_DP
! NAMELIST の読み込み
! NAMELIST is input
!
!!$ if ( trim(namelist_filename) /= '' ) then
!!$ call FileOpen( unit_nml, & ! (out)
!!$ & namelist_filename, mode = 'r' ) ! (in)
!!$
!!$ rewind( unit_nml )
!!$ read( unit_nml, & ! (in)
!!$ & nml = vdiffusion_my_nml, & ! (out)
!!$ & iostat = iostat_nml ) ! (out)
!!$ close( unit_nml )
!!$
!!$ call NmlutilMsg( iostat_nml, module_name ) ! (in)
!!$ if ( iostat_nml == 0 ) write( STDOUT, nml = vdiffusion_my_nml )
!!$ end if
! 補助的な変数を計算するサブルーチン・関数群
! Subroutines and functions for calculating auxiliary variables
!
call AuxVarsInit
#ifdef JMAVDIFF
PressRef = RefPressForPotTemp
TimeStep = 2.0_DP * DelTime
SurfEmis = 1.0_DP
LevelNum = mymodel25
!!$ ntime = nint(itend * timestep / dt_str)
!!$ call pp_monit_ini(.true., timestep, dt_str_in = dt_str)
!!$ call pp_monit_grads_ini(1, 1, pbl_nz, ntime, dt_str, idate, z_f, 'output')
call pp_phys_const_set
call pbl_const_ini( PressRef, TimeStep, SurfEmis )
call pbl_parm_ini
call pbl_grid_ini( kmax )
call pbl_mym_option_ini(levflag_in = LevelNum)
call pbl_mym_parm_ini
call pbl_mym_const_ini
#else
call MessageNotify( 'E', module_name, 'JMA MY library is not included.' )
#endif
! ヒストリデータ出力のためのへの変数登録
! Register of variables for history data output
!
!!$ call HistoryAutoAddVariable( 'VelDiffCoef', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'diffusion coef. momentum', 'm2 s-1' )
!!$ call HistoryAutoAddVariable( 'TempDiffCoef', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'diffusion coef. heat ', 'm2 s-1' )
!!$ call HistoryAutoAddVariable( 'QVapDiffCoef', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'diffusion coef. moisture', 'm2 s-1' )
!!$
!!$ call HistoryAutoAddVariable( 'MomFluxX', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'eastward momentum flux', 'N m-2' )
!!$ call HistoryAutoAddVariable( 'MomFluxY', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'northward momentum flux', 'N m-2' )
!!$ call HistoryAutoAddVariable( 'HeatFlux', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'heat flux', 'W m-2' )
!!$ call HistoryAutoAddVariable( 'QVapFlux', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'moisture flux', 'W m-2' )
!!$
!!$ call HistoryAutoAddVariable( 'DUDtVDiff', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'tendency of zonal wind by vertical diffusion', 'm s-2' )
!!$ call HistoryAutoAddVariable( 'DVDtVDiff', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'tendency of meridional wind by vertical diffusion', 'm s-2' )
!!$ call HistoryAutoAddVariable( 'DTempDtVDiff', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'tendency of temperature by vertical diffusion', 'K s-1' )
!!$ call HistoryAutoAddVariable( 'DQVapDtVDiff', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'tendency of specific humidity by vertical diffusion', 's-1' )
!!$
!!$ call HistoryAutoAddVariable( 'TurKinEne', &
!!$ & (/ AxNameX, AxNameY, AxNameR, AxNameT /), &
!!$ & 'turbulent kinetic energy', 'm2 s-2' )
!!$
!!$ call HistoryAutoAddVariable( 'TKEPShear', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'turbulent kinetic energy production rate by shear', 'm2 s-3' )
!!$ call HistoryAutoAddVariable( 'TKEPBuoy', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'turbulent kinetic energy production rate by buoyancy', 'm2 s-3' )
!!$ call HistoryAutoAddVariable( 'TKEDiss', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'turbulent kinetic energy dissipation rate', 'm2 s-3' )
!!$ call HistoryAutoAddVariable( 'MixLength', &
!!$ & (/ AxNameX, AxNameY, AxNameZ, AxNameT /), &
!!$ & 'mixing length', 'm' )
! 印字 ; Print
!
call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
#ifdef JMAVDIFF
call MessageNotify( 'M', module_name, 'Level MY2.5: %d, MY3: %d', i = (/ mymodel25, mymodel3 /) )
#endif
call MessageNotify( 'M', module_name, ' %d is specified.', i = (/ LevelNum /) )
!!$ call MessageNotify( 'M', module_name, 'For vertical diffusion flux:' )
!!$ call MessageNotify( 'M', module_name, ' FlagConstDiffCoef = %b', l = (/ FlagConstDiffCoef /) )
!!$ call MessageNotify( 'M', module_name, ' ConstDiffCoefM = %f', d = (/ ConstDiffCoefM /) )
!!$ call MessageNotify( 'M', module_name, ' ConstDiffCoefH = %f', d = (/ ConstDiffCoefH /) )
!!$ call MessageNotify( 'M', module_name, ' SquareVelMin = %f', d = (/ SquareVelMin /) )
!!$ call MessageNotify( 'M', module_name, ' BulkRiNumMin = %f', d = (/ BulkRiNumMin /) )
!!$ call MessageNotify( 'M', module_name, 'For diffusion coefficients:' )
!!$ call MessageNotify( 'M', module_name, ' MixLengthMax = %f', d = (/ MixLengthMax /) )
!!$ call MessageNotify( 'M', module_name, ' ShMin = %f', d = (/ ShMin /) )
!!$ call MessageNotify( 'M', module_name, ' SmMin = %f', d = (/ SmMin /) )
!!$ call MessageNotify( 'M', module_name, ' VelDiffCoefMin = %f', d = (/ VelDiffCoefMin /) )
!!$ call MessageNotify( 'M', module_name, ' TempDiffCoefMin = %f', d = (/ TempDiffCoefMin /) )
!!$ call MessageNotify( 'M', module_name, ' VelDiffCoefMax = %f', d = (/ VelDiffCoefMax /) )
!!$ call MessageNotify( 'M', module_name, ' TempDiffCoefMax = %f', d = (/ TempDiffCoefMax /) )
!!$ call MessageNotify( 'M', module_name, ' MYConstA1 = %f', d = (/ MYConstA1 /) )
!!$ call MessageNotify( 'M', module_name, ' MYConstB1 = %f', d = (/ MYConstB1 /) )
!!$ call MessageNotify( 'M', module_name, ' MYConstA2 = %f', d = (/ MYConstA2 /) )
!!$ call MessageNotify( 'M', module_name, ' MYConstB2 = %f', d = (/ MYConstB2 /) )
!!$ call MessageNotify( 'M', module_name, ' MYConstC1 = %f', d = (/ MYConstC1 /) )
call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )
vdiffusion_jma_my_wrap_inited = .true.
end subroutine VDiffusionJMAInit
| Subroutine : | |||
| xyz_U(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyz_V(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax) : | real(DP), intent(in)
| ||
| xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Temp(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_VirTemp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_VirTemp(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_Press(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Press(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xy_SurfHeight(0:imax-1,1:jmax) : | real(DP), intent(in)
| ||
| xyz_Height(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Height(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xyz_Exner(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xyr_Exner(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
| xy_SurfMOLength(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
| xyz_TurKinEne(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
| xy_SurfMomFluxX(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
| xy_SurfMomFluxY(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
| xy_SurfHeatFlux(0:imax-1, 1:jmax) : | real(DP), intent(in) | ||
| xyf_SurfQMixFlux(0:imax-1, 1:jmax, 1:ncmax) : | real(DP), intent(in) | ||
| xyr_MomFluxX(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyr_MomFluxY(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyr_HeatFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyrf_QMixFlux(0:imax-1, 1:jmax, 0:kmax, 1:ncmax) : | real(DP), intent(out)
| ||
| xyr_VelDiffCoef(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyr_TempDiffCoef(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyr_QMixDiffCoef(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
| xyz_DTurKinEneDt(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(out)
|
鉛直拡散フラックスを計算します.
Vertical diffusion flux is calculated by use of MY2.5 model.
subroutine VDiffusionJMAMYWrapper3D( xyz_U, xyz_V, xyzf_QMix, xyz_Temp, xyr_Temp, xyz_VirTemp, xyr_VirTemp, xyz_Press, xyr_Press, xy_SurfHeight, xyz_Height, xyr_Height, xyz_Exner, xyr_Exner, xy_SurfMOLength, xyz_TurKinEne, xy_SurfMomFluxX, xy_SurfMomFluxY, xy_SurfHeatFlux, xyf_SurfQMixFlux, xyr_MomFluxX, xyr_MomFluxY, xyr_HeatFlux, xyrf_QMixFlux, xyr_VelDiffCoef, xyr_TempDiffCoef, xyr_QMixDiffCoef, xyz_DTurKinEneDt )
!
! 鉛直拡散フラックスを計算します.
!
! Vertical diffusion flux is calculated by use of MY2.5 model.
!
! モジュール引用 ; USE statements
!
#ifdef JMAVDIFF
! 物理定数設定
! Physical constants settings
!
use constants, only: GasRDry, CpDry
! $ C_p $ [J kg-1 K-1].
! 乾燥大気の定圧比熱.
! Specific heat of air at constant pressure
#endif
! 時刻管理
! Time control
!
use timeset, only: #ifdef JMAVDIFF #endif use gtool_historyauto, only: HistoryAutoPut
#ifdef JMAVDIFF
!
! JMA library
!
use pbl_grid, only: pbl_grid_set
use pbl_mym_main, only: pbl_mym_main_level25
#endif
! 宣言文 ; Declaration statements
!
implicit none
real(DP), intent(in):: xyz_U (0:imax-1, 1:jmax, 1:kmax)
! $ u $ . 東西風速. Eastward wind
real(DP), intent(in):: xyz_V (0:imax-1, 1:jmax, 1:kmax)
! $ v $ . 南北風速. Northward wind
real(DP), intent(in):: xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax)
! $ q $ . 質量混合比. Mass mixing ratio
real(DP), intent(in):: xyz_Temp (0:imax-1, 1:jmax, 1:kmax)
! $ T $ . 温度. Temperature
real(DP), intent(in):: xyr_Temp (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{T} $ . 温度 (半整数レベル).
! Temperature (half level)
real(DP), intent(in):: xyz_VirTemp (0:imax-1, 1:jmax, 1:kmax)
! $ T_v $ . 仮温度. Virtual temperature
real(DP), intent(in):: xyr_VirTemp (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{T}_v $ . 仮温度 (半整数レベル).
! Virtual temperature (half level)
real(DP), intent(in):: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
! $ p $ . 気圧 (整数レベル).
! Air pressure (full level)
real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
! $ \hat{p} $ . 気圧 (半整数レベル).
! Air pressure (half level)
real(DP), intent(in):: xy_SurfHeight (0:imax-1,1:jmax)
! $ z_s $ . 地表面高度.
! Surface height.
real(DP), intent(in):: xyz_Height (0:imax-1, 1:jmax, 1:kmax)
! 高度 (整数レベル).
! Height (full level)
real(DP), intent(in):: xyr_Height (0:imax-1, 1:jmax, 0:kmax)
! 高度 (半整数レベル).
! Height (half level)
real(DP), intent(in):: xyz_Exner (0:imax-1, 1:jmax, 1:kmax)
! Exner 関数 (整数レベル).
! Exner function (full level)
real(DP), intent(in):: xyr_Exner (0:imax-1, 1:jmax, 0:kmax)
! Exner 関数 (半整数レベル).
! Exner function (half level)
real(DP), intent(in):: xy_SurfMOLength(0:imax-1, 1:jmax)
!
! Turbulent kinetic energy (m2 s-2)
real(DP), intent(in):: xyz_TurKinEne(0:imax-1, 1:jmax, 1:kmax)
!
! Turbulent kinetic energy (m2 s-2)
real(DP), intent(in):: xy_SurfMomFluxX (0:imax-1, 1:jmax)
!
! Eastward momentum flux at surface
real(DP), intent(in):: xy_SurfMomFluxY (0:imax-1, 1:jmax)
!
! Northward momentum flux at surface
real(DP), intent(in):: xy_SurfHeatFlux (0:imax-1, 1:jmax)
real(DP), intent(in):: xyf_SurfQMixFlux(0:imax-1, 1:jmax, 1:ncmax)
real(DP), intent(out):: xyr_MomFluxX (0:imax-1, 1:jmax, 0:kmax)
! 東西方向運動量フラックス.
! Eastward momentum flux
real(DP), intent(out):: xyr_MomFluxY (0:imax-1, 1:jmax, 0:kmax)
! 南北方向運動量フラックス.
! Northward momentum flux
real(DP), intent(out):: xyr_HeatFlux (0:imax-1, 1:jmax, 0:kmax)
! 熱フラックス.
! Heat flux
real(DP), intent(out):: xyrf_QMixFlux(0:imax-1, 1:jmax, 0:kmax, 1:ncmax)
! 質量フラックス.
! Mass flux of compositions
real(DP), intent(out):: xyr_VelDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:運動量.
! Diffusion coefficient: velocity
real(DP), intent(out):: xyr_TempDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:温度.
! Diffusion coefficient: temperature
real(DP), intent(out):: xyr_QMixDiffCoef (0:imax-1, 1:jmax, 0:kmax)
! 拡散係数:比湿.
! Diffusion coefficient: specific humidity
real(DP), intent(out):: xyz_DTurKinEneDt (0:imax-1, 1:jmax, 1:kmax)
!
! Tendency of turbulent kinetic energy
! 作業変数
! Work variables
!
#ifdef JMAVDIFF
real(DP) :: z_U (1:kmax)
real(DP) :: z_V (1:kmax)
real(DP) :: z_QH2OVap(1:kmax)
real(DP) :: z_QH2OLiq(1:kmax)
real(DP) :: z_QH2OSol(1:kmax)
real(DP) :: z_Temp (1:kmax)
real(DP) :: r_Temp (0:kmax)
real(DP) :: z_VirTemp (1:kmax)
real(DP) :: r_VirTemp (0:kmax)
real(DP) :: z_Press (1:kmax)
real(DP) :: r_Press (0:kmax)
real(DP) :: SurfHeight
real(DP) :: z_Height (1:kmax)
real(DP) :: r_Height (0:kmax)
real(DP) :: z_Exner (1:kmax)
real(DP) :: r_Exner (0:kmax)
real(DP) :: z_PotTemp(1:kmax)
real(DP) :: z_QCldWat(1:kmax)
real(DP) :: z_QCldIce(1:kmax)
real(DP) :: z_TurKinEne(1:kmax)
real(DP) :: SurfMomFluxX
real(DP) :: SurfMomFluxY
real(DP) :: SurfMomFlux
real(DP) :: SurfHeatFlux
real(DP) :: SurfQVapFlux
real(DP) :: r_MomFluxX (0:kmax)
real(DP) :: r_MomFluxY (0:kmax)
real(DP) :: r_HeatFlux (0:kmax)
real(DP) :: rf_QMixFlux(0:kmax, 1:ncmax)
real(DP) :: r_VelDiffCoef (0:kmax)
real(DP) :: r_TempDiffCoef(0:kmax)
real(DP) :: r_QMixDiffCoef(0:kmax)
real(DP) :: z_DTurKinEneDt (1:kmax)
!
! Tendency of turbulent kinetic energy
real(DP) :: r_Rho(0:kmax)
real(DP) :: SurfMOLengthInv
real(DP) :: z_DelHeight (1:kmax-1)
real(DP) :: z_DelHeightInv(1:kmax)
real(DP) :: r_DelHeightInv(1:kmax-1)
real(DP) :: r_LowWeightFactF2H(1:kmax-1)
real(DP) :: r_UppWeightFactF2H(1:kmax-1)
real(DP) :: r_LowWeightFactH2F(1:kmax-1)
real(DP) :: r_UppWeightFactH2F(1:kmax-1)
real(DP) :: z_QKE (1:kmax)
real(DP) :: z_TSq (1:kmax)
real(DP) :: z_QSq (1:kmax)
real(DP) :: z_Cov (1:kmax)
real(DP) :: z_TSqA(1:kmax)
real(DP) :: z_QSqA(1:kmax)
real(DP) :: z_CovA(1:kmax)
integer:: i ! 経度方向に回る DO ループ用作業変数
! Work variables for DO loop in longitude
integer:: j ! 緯度方向に回る DO ループ用作業変数
! Work variables for DO loop in latitude
integer:: n
#endif
! 実行文 ; Executable statement
!
! 初期化確認
! Initialization check
!
if ( .not. vdiffusion_jma_my_wrap_inited ) then
call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
end if
#ifdef JMAVDIFF
! 計算時間計測開始
! Start measurement of computation time
!
call TimesetClockStart( module_name )
do j = 1, jmax
do i = 0, imax-1
z_U = xyz_U (i,j,:)
z_V = xyz_V (i,j,:)
z_QH2OVap = xyzf_QMix (i,j,:,IndexH2OVap)
if ( IndexH2OLiq > 0 ) then
z_QCldWat = xyzf_QMix (i,j,:,IndexH2OLiq)
else
z_QCldWat = 0.0_DP
end if
if ( IndexH2OSol > 0 ) then
z_QCldIce = xyzf_QMix (i,j,:,IndexH2OSol)
else
z_QCldIce = 0.0_DP
end if
z_Temp = xyz_Temp (i,j,:)
r_Temp = xyr_Temp (i,j,:)
z_VirTemp = xyz_VirTemp (i,j,:)
r_VirTemp = xyr_VirTemp (i,j,:)
z_Press = xyz_Press (i,j,:)
r_Press = xyr_Press (i,j,:)
SurfHeight = xy_SurfHeight (i,j)
z_Height = xyz_Height (i,j,:) - xyr_Height (i,j,0)
r_Height = xyr_Height (i,j,:) - xyr_Height (i,j,0)
z_Exner = xyz_Exner (i,j,:)
r_Exner = xyr_Exner (i,j,:)
z_TurKinEne = xyz_TurKinEne (i,j,:)
r_Rho = r_Press / ( GasRDry * r_VirTemp )
SurfMomFluxX = xy_SurfMomFluxX(i,j) / r_Rho(0)
SurfMomFluxY = xy_SurfMomFluxY(i,j) / r_Rho(0)
SurfHeatFlux = xy_SurfHeatFlux(i,j) / ( r_Rho(0) * CpDry )
SurfQVapFlux = xyf_SurfQMixFlux(i,j,IndexH2OVap) / r_Rho(0)
SurfMOLengthInv = 1.0_DP / xy_SurfMOLength(i,j)
z_PotTemp = z_Temp / z_Exner
SurfMomFlux = sqrt( SurfMomFluxX**2 + SurfMomFluxY**2 )
call pbl_grid_set( z_Height, r_Height(1:kmax-1), z_DelHeight(1:kmax-1), z_DelHeightInv(1:kmax), r_DelHeightInv(1:kmax-1), r_LowWeightFactF2H(1:kmax-1), r_UppWeightFactF2H(1:kmax-1), r_LowWeightFactH2F(1:kmax-1), r_UppWeightFactH2F(1:kmax-1) )
!!$ if ( FlagFirst ) then
!!$ call pbl_mym_initialize_run(&
!!$ & ftl_surf_ex, fqw_surf_ex, uf, l_mo_inv, &
!!$ & u, v, pt, qv, qc, qci, prs, &
!!$ & z_f, dz_f, rdz_f, rdz_h, f2h_m, f2h_p, &
!!$ & qke, tsq, qsq, cov)
!!$ call pbl_mym_initialize_run(&
!!$ & r_HeatFlux(0), rf_QMixFlux(0,IndexH2OVap), & ! (in )
!!$ & SurfMomFlux, l_mo_inv, & ! (in )
!!$ & z_U, z_V, z_PotTemp, zf_QMix(:,IndexH2OVap), & ! (in )
!!$ & z_QCldWat, z_QCldIce, z_Press, & ! (in )
!!$ & z_Height, z_DelHeight(1:kmax-1), & ! (in )
!!$ & z_DelHeightInv(1:kmax), r_DelHeightInv(1:kmax-1), & ! (in )
!!$ & r_LowWeightFactF2H(1:kmax-1), r_UppWeightFactF2H(1:kmax-1), & ! (in )
!!$ & qke, tsq, qsq, cov & ! (out)
!!$ & )
!!$ FlagFirst = .false.
!!$ end if
z_QKE = z_TurKinEne * 2.0_DP
z_TSq = 0.0_DP
z_QSq = 0.0_DP
z_Cov = 0.0_DP
call pbl_mym_main_level25( SurfMomFluxX, SurfMomFluxY, SurfHeatFlux, SurfQVapFlux, SurfMomFlux, SurfMOLengthInv, r_Press(0), z_U, z_V, z_PotTemp, z_QH2OVap, z_QCldWat, z_QCldIce, z_Press, z_Exner, z_QKE, z_TSq, z_QSq, z_Cov, z_Height, z_DelHeight(1:kmax-1), z_DelHeightInv(1:kmax), r_DelHeightInv(1:kmax-1), r_LowWeightFactF2H(1:kmax-1), r_UppWeightFactF2H(1:kmax-1), r_LowWeightFactH2F(1:kmax-1), r_UppWeightFactH2F(1:kmax-1), r_VelDiffCoef(1:kmax), r_TempDiffCoef(1:kmax), z_DTurKinEneDt, z_TSqA, z_QSqA, z_CovA, r_MomFluxX(1:kmax), r_MomFluxY(1:kmax), r_HeatFlux(1:kmax), rf_QMixFlux(1:kmax,IndexH2OVap) )
r_VelDiffCoef (0) = 0.0_DP
r_TempDiffCoef(0) = 0.0_DP
z_DTurKinEneDt = z_DTurKinEneDt / 2.0_DP
r_MomFluxX(0) = 0.0_DP
r_MomFluxY(0) = 0.0_DP
r_HeatFlux(0) = 0.0_DP
rf_QMixFlux(0,IndexH2OVap) = 0.0_DP
r_QMixDiffCoef = r_TempDiffCoef
do n = 1, IndexH2OVap-1
rf_QMixFlux(:,n) = 0.0_DP
end do
do n = IndexH2OVap+1, ncmax
rf_QMixFlux(:,n) = 0.0_DP
end do
r_MomFluxX = r_Rho * r_MomFluxX
r_MomFluxY = r_Rho * r_MomFluxY
r_HeatFlux = r_Rho * CpDry * r_HeatFlux * r_Exner
rf_QMixFlux(:,IndexH2OVap) = r_Rho * rf_QMixFlux(:,IndexH2OVap)
!!$ call VDiffusionMY251D( &
!!$ & z_U, z_V, zf_QMix, & ! (in)
!!$ & z_Temp, r_Temp, z_VirTemp, r_VirTemp, r_Press, & ! (in)
!!$ & SurfHeight, & ! (in)
!!$ & z_Height, r_Height, z_Exner, r_Exner, & ! (in)
!!$ & z_TurKinEne, & ! (in)
!!$ & SurfMomFluxX, SurfMomFluxY, & ! (in)
!!$ & r_MomFluxX, r_MomFluxY, r_HeatFlux, rf_QMixFlux, & ! (out)
!!$ & r_VelDiffCoef, r_TempDiffCoef, r_QMixDiffCoef, & ! (out)
!!$ & z_DTurKinEneDt & ! (out)
!!$ & )
xyr_MomFluxX (i,j,:) = r_MomFluxX
xyr_MomFluxY (i,j,:) = r_MomFluxY
xyr_HeatFlux (i,j,:) = r_HeatFlux
xyrf_QMixFlux (i,j,:,:) = rf_QMixFlux
xyr_VelDiffCoef (i,j,:) = r_VelDiffCoef
xyr_TempDiffCoef(i,j,:) = r_TempDiffCoef
xyr_QMixDiffCoef(i,j,:) = r_QMixDiffCoef
xyz_DTurKinEneDt(i,j,:) = z_DTurKinEneDt
end do
end do
#else
call MessageNotify( 'E', module_name, 'JMA MY library is not included.' )
#endif
! 拡散係数の出力
! Output diffusion coefficients
!
! 拡散係数出力
! Diffusion coeffficients output
!
!!$ call HistoryAutoPut( TimeN, 'VelDiffCoef', xyr_VelDiffCoef )
!!$ call HistoryAutoPut( TimeN, 'TempDiffCoef', xyr_TempDiffCoef )
!!$ call HistoryAutoPut( TimeN, 'QVapDiffCoef', xyr_QMixDiffCoef )
!!$
!!$ call HistoryAutoPut( TimeN, 'TKEPShear', xyz_TurKinEneProShear )
!!$ call HistoryAutoPut( TimeN, 'TKEPBuoy' , xyz_TurKinEneProBuoya )
!!$ call HistoryAutoPut( TimeN, 'TKEDiss' , xyz_TurKinEneDiss )
!!$
!!$ call HistoryAutoPut( TimeN, 'MixLength' , xyz_MixLength )
! 計算時間計測一時停止
! Pause measurement of computation time
!
call TimesetClockStop( module_name )
end subroutine VDiffusionJMAMYWrapper3D
| Constant : | |||
| module_name = ‘vdiffusion_jma_my_wrapper‘ : | character(*), parameter
|
| Variable : | |||
| vdiffusion_jma_my_wrap_inited = .false. : | logical, save
|
| Constant : | |||
| version = ’$Name: $’ // ’$Id: vdiffusion_jma_my_wrapper.F90,v 1.2 2015/02/11 11:53:55 yot Exp $’ : | character(*), parameter
|