Class | rad_rte_nonscat |
In: |
radiation/rad_rte_nonscat.f90
|
Note that Japanese and English are described in parallel.
!$ ! RadDTempDt : | 放射フラックスによる温度変化の計算 |
!$ ! RadFluxOutput : | 放射フラックスの出力 |
!$ ! ———— : | ———— |
!$ ! RadDTempDt : | Calculate temperature tendency with radiation flux |
!$ ! RadFluxOutput : | Output radiation fluxes |
Subroutine : | |||
xyz_IntPF(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
xy_SurfIntPF(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_IntDPFDT1(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_SurfIntDPFDT(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xyrr_Trans(0:imax-1, 1:jmax, 0:kmax, 0:kmax) : | real(DP), intent(in )
| ||
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
|
散乱なしの場合の放射伝達方程式の計算
Integrate radiative transfer equation without scattering
subroutine RadRTENonScat( xyz_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_SurfIntDPFDT, xyrr_Trans, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux ) ! ! 散乱なしの場合の放射伝達方程式の計算 ! ! Integrate radiative transfer equation without scattering ! ! モジュール引用 ; USE statements ! ! 宣言文 ; Declaration statements ! real(DP), intent(in ) :: xyz_IntPF (0:imax-1, 1:jmax, 1:kmax) ! Integrated Planck function real(DP), intent(in ) :: xy_SurfIntPF (0:imax-1, 1:jmax) ! Integrated Planck function with surface temperature real(DP), intent(in ) :: xy_IntDPFDT1 (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function real(DP), intent(in ) :: xy_SurfIntDPFDT (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function ! with surface temperature real(DP), intent(in ) :: xyrr_Trans (0:imax-1, 1:jmax, 0:kmax, 0:kmax) ! 透過係数. ! Transmission coefficient real(DP), intent(out) :: xyr_RadLUwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Upward longwave flux real(DP), intent(out) :: xyr_RadLDwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Downward longwave flux real(DP), intent(out) :: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Upward longwave flux derivative real(DP), intent(out) :: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Downward longwave flux derivative ! 作業変数 ! Work variables ! integer:: k, kk ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. rad_rte_nonscat_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! 放射フラックス計算 ! Calculate radiation flux ! ! Initialization ! xyr_RadLDwFlux = 0.0_DP xyr_RadLUwFlux = 0.0_DP ! ! Downward flux ! do k = kmax, 0, -1 do kk = kmax, k+1, -1 xyr_RadLDwFlux(:,:,k) = xyr_RadLDwFlux(:,:,k) + xyz_IntPF(:,:,kk) * ( xyrr_Trans(:,:,k,kk-1) - xyrr_Trans(:,:,k,kk) ) end do end do ! ! Upward flux ! ! Set upward flux ! do k = 0, kmax xyr_RadLUwFlux(:,:,k) = xy_SurfIntPF * xyrr_Trans(:,:,k,0) do kk = 1, k xyr_RadLUwFlux(:,:,k) = xyr_RadLUwFlux(:,:,k) - xyz_IntPF(:,:,kk) * ( xyrr_Trans(:,:,k,kk-1) - xyrr_Trans(:,:,k,kk) ) end do end do ! 放射フラックスの変化率の計算 ! Calculate rate of change of radiative flux ! do k = 0, kmax xyra_DelRadLUwFlux(:,:,k,0) = xy_SurfIntDPFDT * xyrr_Trans(:,:,k,0) end do k = 0 xyra_DelRadLUwFlux(:,:,k,1) = 0.0_DP do k = 1, kmax xyra_DelRadLUwFlux(:,:,k,1) = - xy_IntDPFDT1(:,:) * ( xyrr_Trans(:,:,k,0) - xyrr_Trans(:,:,k,1) ) end do do k = 0, kmax xyra_DelRadLDwFlux(:,:,k,0) = 0.0_DP end do k = 0 xyra_DelRadLDwFlux(:,:,k,1) = + xy_IntDPFDT1(:,:) * ( xyrr_Trans(:,:,k,0) - xyrr_Trans(:,:,k,1) ) do k = 1, kmax xyra_DelRadLDwFlux(:,:,k,1) = 0.0_DP end do end subroutine RadRTENonScat
Subroutine : |
rad_rte_nonscat モジュールの初期化を行います. NAMELIST#rad_rte_nonscat_nml の読み込みはこの手続きで行われます.
"rad_rte_nonscat" module is initialized. "NAMELIST#rad_rte_nonscat_nml" is loaded in this procedure.
This procedure input/output NAMELIST#rad_rte_nonscat_nml .
subroutine RadRTENonScatInit ! ! rad_rte_nonscat モジュールの初期化を行います. ! NAMELIST#rad_rte_nonscat_nml の読み込みはこの手続きで行われます. ! ! "rad_rte_nonscat" module is initialized. ! "NAMELIST#rad_rte_nonscat_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 ! ヒストリデータ出力 ! History data output ! use gtool_historyauto, only: HistoryAutoAddVariable ! 文字列操作 ! Character handling ! use dc_string, only: toChar ! ガウス重み, 分点の計算 ! Calculate Gauss node and Gaussian weight ! use gauss_quad, only : GauLeg ! 宣言文 ; Declaration statements ! integer:: unit_nml ! NAMELIST ファイルオープン用装置番号. ! Unit number for NAMELIST file open integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT. ! IOSTAT of NAMELIST read ! NAMELIST 変数群 ! NAMELIST group name ! namelist /rad_rte_nonscat_nml/ DiffFact, NumGaussNodeZAInt ! ! デフォルト値については初期化手続 "rad_rte_nonscat#RadRTENonScatInit" ! のソースコードを参照のこと. ! ! Refer to source codes in the initialization procedure ! "rad_rte_nonscat#RadRTENonScatInit" for the default values. ! ! 実行文 ; Executable statement ! if ( rad_rte_nonscat_inited ) return ! デフォルト値の設定 ! Default values settings ! NumGaussNodeZAInt = 3 DiffFact = 1.66_DP ! NAMELIST の読み込み ! NAMELIST is input ! if ( trim(namelist_filename) /= '' ) then call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in) rewind( unit_nml ) read( unit_nml, nml = rad_rte_nonscat_nml, iostat = iostat_nml ) ! (out) close( unit_nml ) call NmlutilMsg( iostat_nml, module_name ) ! (in) end if if ( NumGaussNodeZAInt > 0 ) then allocate( a_CosZA ( NumGaussNodeZAInt ) ) allocate( a_GaussWeight( NumGaussNodeZAInt ) ) call GauLeg( 0.0_DP, 1.0_DP, NumGaussNodeZAInt, a_CosZA, a_GaussWeight ) end if ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, 'NumGaussNodeZAInt = %d', i = (/ NumGaussNodeZAInt /) ) call MessageNotify( 'M', module_name, 'DiffFact = %f', d = (/ DiffFact /) ) call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) ) rad_rte_nonscat_inited = .true. end subroutine RadRTENonScatInit
Subroutine : | |||
xy_SurfAlbedo(0:imax-1, 1:jmax) : | real(DP), intent(in ) | ||
xyr_IntPF(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in )
| ||
xy_SurfIntPF(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_IntDPFDT1(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_SurfIntDPFDT(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xyz_OptDep(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
|
散乱なしの場合の放射伝達方程式の計算
Integrate radiative transfer equation without scattering
subroutine RadRTENonScatMonoSemiAnal( xy_SurfAlbedo, xyr_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_SurfIntDPFDT, xyz_OptDep, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux ) ! ! 散乱なしの場合の放射伝達方程式の計算 ! ! Integrate radiative transfer equation without scattering ! ! モジュール引用 ; USE statements ! ! 宣言文 ; Declaration statements ! real(DP), intent(in ) :: xy_SurfAlbedo (0:imax-1, 1:jmax) real(DP), intent(in ) :: xyr_IntPF (0:imax-1, 1:jmax, 0:kmax) ! Integrated Planck function real(DP), intent(in ) :: xy_SurfIntPF (0:imax-1, 1:jmax) ! Integrated Planck function with surface temperature real(DP), intent(in ) :: xy_IntDPFDT1 (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function real(DP), intent(in ) :: xy_SurfIntDPFDT (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function ! with surface temperature real(DP), intent(in ) :: xyz_OptDep (0:imax-1, 1:jmax, 1:kmax) ! 光学的厚さ. ! Optical depth real(DP), intent(out) :: xyr_RadLUwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Upward longwave flux real(DP), intent(out) :: xyr_RadLDwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Downward longwave flux real(DP), intent(out) :: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Upward longwave flux derivative real(DP), intent(out) :: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Downward longwave flux derivative ! 作業変数 ! Work variables ! real(DP) :: CosZA real(DP) :: GaussWeight real(DP) :: xyz_TransEachLayer0(0:imax-1, 1:jmax, 1:kmax) real(DP) :: xyz_TransEachLayer (0:imax-1, 1:jmax, 1:kmax) ! 透過係数. ! Transmission coefficient real(DP) :: xyr_Trans0(0:imax-1, 1:jmax, 0:kmax) ! ! Transmission coefficient from surface to layer interfaces real(DP) :: xyr_Trans1(0:imax-1, 1:jmax, 0:kmax) ! ! Transmission coefficient from top of the lowest layer to layer interfaces real(DP) :: xyz_DPFDOptDep(0:imax-1, 1:jmax, 1:kmax) real(DP) :: xyr_UwFlux (0:imax-1, 1:jmax, 0:kmax) real(DP) :: xyr_DwFlux (0:imax-1, 1:jmax, 0:kmax) real(DP) :: xyra_DelUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) real(DP) :: xyra_DelDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) real(DP) :: IntFact integer:: k ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction integer:: n ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. rad_rte_nonscat_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! 放射フラックス計算 ! Calculate radiation flux ! xyz_TransEachLayer0 = exp( - xyz_OptDep ) do k = 1, kmax xyz_DPFDOptDep(:,:,k) = ( xyr_IntPF(:,:,k-1) - xyr_IntPF(:,:,k) ) / max( xyz_OptDep(:,:,k), 1.0d-100 ) end do if ( NumGaussNodeZAInt > 0 ) then ! for case with Gaussian quadrature IntFact = 2.0_DP else ! for two-stream approximation or case with diffusivity factor IntFact = 1.0_DP end if ! Initialization ! xyr_RadLUwFlux = 0.0_DP xyr_RadLDwFlux = 0.0_DP xyra_DelRadLUwFlux = 0.0_DP xyra_DelRadLDwFlux = 0.0_DP ! Loop for Gaussian quadrature ! loop_gq : do n = 1, max( NumGaussNodeZAInt, 1 ) ! Preparetion ! if ( NumGaussNodeZAInt > 0 ) then CosZA = a_CosZA(n) else CosZA = 1.0_DP / DiffFact end if ! xyz_TransEachLayer = ( xyz_TransEachLayer0 )**(1.0d0/CosZA) ! xyr_Trans0(:,:,0) = 1.0_DP do k = 1, kmax xyr_Trans0(:,:,k) = xyr_Trans0(:,:,k-1) * xyz_TransEachLayer(:,:,k) end do xyr_Trans1(:,:,0) = xyz_TransEachLayer(:,:,1) xyr_Trans1(:,:,1) = 1.0_DP do k = 2, kmax xyr_Trans1(:,:,k) = xyr_Trans0(:,:,k-1) * xyz_TransEachLayer(:,:,k) end do ! ! Downward flux ! k = kmax xyr_DwFlux(:,:,k) = 0.0_DP do k = kmax-1, 0, -1 xyr_DwFlux(:,:,k) = ( xyr_DwFlux(:,:,k+1) - IntFact * xyr_IntPF(:,:,k+1) ) * xyz_TransEachLayer(:,:,k+1) + IntFact * xyr_IntPF(:,:,k) - IntFact * CosZA * xyz_DPFDOptDep(:,:,k+1) * ( 1.0_DP - xyz_TransEachLayer(:,:,k+1) ) end do ! ! Upward flux ! k = 0 xyr_UwFlux(:,:,k) = IntFact * xy_SurfIntPF + xy_SurfAlbedo * xyr_DwFlux(:,:,0) do k = 1, kmax xyr_UwFlux(:,:,k) = ( xyr_UwFlux(:,:,k-1) - IntFact * xyr_IntPF(:,:,k-1) ) * xyz_TransEachLayer(:,:,k) + IntFact * xyr_IntPF(:,:,k) + IntFact * CosZA * xyz_DPFDOptDep(:,:,k) * ( 1.0_DP - xyz_TransEachLayer(:,:,k) ) end do ! 放射フラックスの変化率の計算 ! Calculate rate of change of radiative flux ! do k = 0, kmax xyra_DelUwFlux(:,:,k,0) = IntFact * xy_SurfIntDPFDT * xyr_Trans0(:,:,k) end do do k = 0, kmax xyra_DelDwFlux(:,:,k,0) = 0.0_DP end do xyra_DelUwFlux(:,:,:,1) = 0.0_DP xyra_DelDwFlux(:,:,:,1) = 0.0_DP ! Sum over zenith angle ! if ( NumGaussNodeZAInt > 0 ) then GaussWeight = a_GaussWeight( n ) xyr_RadLUwFlux = xyr_RadLUwFlux + xyr_UwFlux * CosZA * GaussWeight xyr_RadLDwFlux = xyr_RadLDwFlux + xyr_DwFlux * CosZA * GaussWeight xyra_DelRadLUwFlux = xyra_DelRadLUwFlux + xyra_DelUwFlux * CosZA * GaussWeight xyra_DelRadLDwFlux = xyra_DelRadLDwFlux + xyra_DelDwFlux * CosZA * GaussWeight else xyr_RadLUwFlux = xyr_UwFlux xyr_RadLDwFlux = xyr_DwFlux xyra_DelRadLUwFlux = xyra_DelUwFlux xyra_DelRadLDwFlux = xyra_DelDwFlux end if end do loop_gq end subroutine RadRTENonScatMonoSemiAnal
Subroutine : | |||
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
xy_SurfTemp(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_SurfEmis(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xyrr_Trans(0:imax-1, 1:jmax, 0:kmax, 0:kmax) : | real(DP), intent(in )
| ||
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
WNs : | real(DP), intent(in ), optional | ||
WNe : | real(DP), intent(in ), optional | ||
NumGaussNode : | integer , intent(in ), optional |
散乱なしの場合の放射伝達方程式の計算
Integrate radiative transfer equation without scattering
subroutine RadRTENonScatWrapper( xyz_Temp, xy_SurfTemp, xy_SurfEmis, xyrr_Trans, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux, WNs, WNe, NumGaussNode ) ! ! 散乱なしの場合の放射伝達方程式の計算 ! ! Integrate radiative transfer equation without scattering ! ! モジュール引用 ; USE statements ! ! プランク関数の計算 ! Calculate Planck function ! use planck_func, only : Integ_PF_GQ_Array3D, Integ_PF_GQ_Array2D, Integ_DPFDT_GQ_Array2D ! 宣言文 ; Declaration statements ! real(DP), intent(in ) :: xyz_Temp (0:imax-1, 1:jmax, 1:kmax) ! $ T $ . 温度. Temperature real(DP), intent(in ) :: xy_SurfTemp (0:imax-1, 1:jmax) ! 地表面温度. ! Surface temperature real(DP), intent(in ) :: xy_SurfEmis (0:imax-1, 1:jmax) ! 惑星表面射出率. ! Surface emissivity real(DP), intent(in ) :: xyrr_Trans (0:imax-1, 1:jmax, 0:kmax, 0:kmax) ! 透過係数. ! Transmission coefficient real(DP), intent(out) :: xyr_RadLUwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Upward longwave flux real(DP), intent(out) :: xyr_RadLDwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Downward longwave flux real(DP), intent(out) :: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Upward longwave flux derivative real(DP), intent(out) :: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Downward longwave flux derivative real(DP), intent(in ), optional :: WNs real(DP), intent(in ), optional :: WNe integer , intent(in ), optional :: NumGaussNode ! 作業変数 ! Work variables ! real(DP):: xyz_IntPF (0:imax-1, 1:jmax, 1:kmax) ! Integrated Planck function real(DP):: xy_SurfIntPF (0:imax-1, 1:jmax) ! Integrated Planck function with surface temperature real(DP):: xy_IntDPFDT1 (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function real(DP):: xy_SurfIntDPFDT (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function ! with surface temperature ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. rad_rte_nonscat_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! Check arguments ! if ( present( WNs ) .or. present( WNe ) .or. present( NumGaussNode ) ) then if ( .not. ( present( WNs ) .and. present( WNe ) .and. present( NumGaussNode ) ) ) then call MessageNotify( 'E', module_name, 'All of WNs, WNe, and NumGaussNode have to be present.' ) end if end if if ( present( WNs ) ) then ! Case for non-grey atmosphere ! ! Integrate Planck function and temperature derivative of it ! call Integ_PF_GQ_Array3D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, 1, kmax, xyz_Temp, xyz_IntPF ) call Integ_PF_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_SurfIntPF ) call Integ_DPFDT_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xyz_Temp(:,:,1), xy_IntDPFDT1 ) call Integ_DPFDT_GQ_Array2D( WNs, WNe, NumGaussNode, 0, imax-1, 1, jmax, xy_SurfTemp, xy_SurfIntDPFDT ) xyz_IntPF = PI * xyz_IntPF xy_SurfIntPF = xy_SurfEmis * PI * xy_SurfIntPF xy_IntDPFDT1 = PI * xy_IntDPFDT1 xy_SurfIntDPFDT = xy_SurfEmis * PI * xy_SurfIntDPFDT else ! Case for grey atmosphere ! xyz_IntPF = StB * xyz_Temp**4 xy_SurfIntPF = xy_SurfEmis * StB * xy_SurfTemp**4 xy_IntDPFDT1 = 4.0_DP * StB * xyz_Temp(:,:,1)**3 xy_SurfIntDPFDT = xy_SurfEmis * 4.0_DP * StB * xy_SurfTemp**3 end if call RadRTENonScat( xyz_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_SurfIntDPFDT, xyrr_Trans, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux ) end subroutine RadRTENonScatWrapper
Variable : | |||
rad_rte_nonscat_inited = .false. : | logical, save, public
|
Subroutine : | |||
DiffFact : | real(DP), intent(in )
| ||
xyz_IntPF(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
xy_SurfIntPF(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_IntDPFDT1(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xy_SurfIntDPFDT(0:imax-1, 1:jmax) : | real(DP), intent(in )
| ||
xyz_OptDep(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in )
| ||
xyr_RadLUwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_RadLDwFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLUwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
xyra_DelRadLDwFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
|
散乱なしの場合の放射伝達方程式の計算
Integrate radiative transfer equation without scattering
subroutine RadRTENonScatMonoWithDiffFact( DiffFact, xyz_IntPF, xy_SurfIntPF, xy_IntDPFDT1, xy_SurfIntDPFDT, xyz_OptDep, xyr_RadLUwFlux, xyr_RadLDwFlux, xyra_DelRadLUwFlux, xyra_DelRadLDwFlux ) ! ! 散乱なしの場合の放射伝達方程式の計算 ! ! Integrate radiative transfer equation without scattering ! ! モジュール引用 ; USE statements ! ! 宣言文 ; Declaration statements ! real(DP), intent(in ) :: DiffFact ! Diffusivity factor real(DP), intent(in ) :: xyz_IntPF (0:imax-1, 1:jmax, 1:kmax) ! Integrated Planck function real(DP), intent(in ) :: xy_SurfIntPF (0:imax-1, 1:jmax) ! Integrated Planck function with surface temperature real(DP), intent(in ) :: xy_IntDPFDT1 (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function real(DP), intent(in ) :: xy_SurfIntDPFDT (0:imax-1, 1:jmax) ! Integrated temperature derivative of Planck function ! with surface temperature real(DP), intent(in ) :: xyz_OptDep (0:imax-1, 1:jmax, 1:kmax) ! 光学的厚さ. ! Optical depth real(DP), intent(out) :: xyr_RadLUwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Upward longwave flux real(DP), intent(out) :: xyr_RadLDwFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Downward longwave flux real(DP), intent(out) :: xyra_DelRadLUwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Upward longwave flux derivative real(DP), intent(out) :: xyra_DelRadLDwFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Downward longwave flux derivative ! 作業変数 ! Work variables ! real(DP) :: xyz_TransEachLayer(0:imax-1, 1:jmax, 1:kmax) ! 透過係数. ! Transmission coefficient real(DP) :: xyr_Trans0(0:imax-1, 1:jmax, 0:kmax) ! ! Transmission coefficient from surface to layer interfaces real(DP) :: xyr_Trans1(0:imax-1, 1:jmax, 0:kmax) ! ! Transmission coefficient from top of the lowest layer to layer interfaces integer:: k ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction ! 実行文 ; Executable statement ! ! 初期化確認 ! Initialization check ! if ( .not. rad_rte_nonscat_inited ) then call MessageNotify( 'E', module_name, 'This module has not been initialized.' ) end if ! 放射フラックス計算 ! Calculate radiation flux ! xyz_TransEachLayer = exp( - DiffFact * xyz_OptDep ) ! xyr_Trans0(:,:,0) = 1.0_DP do k = 1, kmax xyr_Trans0(:,:,k) = xyr_Trans0(:,:,k-1) * xyz_TransEachLayer(:,:,k) end do xyr_Trans1(:,:,0) = xyz_TransEachLayer(:,:,1) xyr_Trans1(:,:,1) = 1.0_DP do k = 2, kmax xyr_Trans1(:,:,k) = xyr_Trans0(:,:,k-1) * xyz_TransEachLayer(:,:,k) end do ! Initialization ! xyr_RadLDwFlux = 0.0_DP xyr_RadLUwFlux = 0.0_DP ! ! Downward flux ! k = kmax xyr_RadLDwFlux(:,:,k) = 0.0_DP do k = kmax-1, 0, -1 xyr_RadLDwFlux(:,:,k) = xyr_RadLDwFlux(:,:,k+1) * xyz_TransEachLayer(:,:,k+1) + xyz_IntPF(:,:,k+1) * ( 1.0_DP - xyz_TransEachLayer(:,:,k+1) ) end do ! ! Upward flux ! k = 0 xyr_RadLUwFlux(:,:,k) = xy_SurfIntPF do k = 1, kmax xyr_RadLUwFlux(:,:,k) = xyr_RadLUwFlux(:,:,k-1) * xyz_TransEachLayer(:,:,k) - xyz_IntPF(:,:,k) * ( xyz_TransEachLayer(:,:,k) - 1.0_DP ) end do ! 放射フラックスの変化率の計算 ! Calculate rate of change of radiative flux ! do k = 0, kmax xyra_DelRadLUwFlux(:,:,k,0) = xy_SurfIntDPFDT * xyr_Trans0(:,:,k) end do k = 0 xyra_DelRadLUwFlux(:,:,k,1) = 0.0_DP do k = 1, kmax xyra_DelRadLUwFlux(:,:,k,1) = - xy_IntDPFDT1 * ( xyr_Trans0(:,:,k) - xyr_Trans1(:,:,k) ) end do do k = 0, kmax xyra_DelRadLDwFlux(:,:,k,0) = 0.0_DP end do k = 0 xyra_DelRadLDwFlux(:,:,k,1) = + xy_IntDPFDT1 * ( xyr_Trans0(:,:,k) - xyr_Trans1(:,:,k) ) do k = 1, kmax xyra_DelRadLDwFlux(:,:,k,1) = 0.0_DP end do end subroutine RadRTENonScatMonoWithDiffFact
Constant : | |||
module_name = ‘rad_rte_nonscat‘ : | character(*), parameter
|
Constant : | |||
version = ’$Name: dcpam5-20130921 $’ // ’$Id: rad_rte_nonscat.f90,v 1.4 2013-05-25 06:49:00 yot Exp $’ : | character(*), parameter
|