Class | rad_DennouAGCM |
In: |
radiation/rad_DennouAGCM.f90
|
Note that Japanese and English are described in parallel.
温度, 比湿, 気圧から, 放射フラックスを計算する放射モデルです.
This is a radiation model that calculates radiation flux from temperature, specific humidity, and air pressure.
Numaguti, Atusi, 1992: 熱帯における積雲活動の大規模構造に関する数値実験. Doctor thesis, 218pp. (Japanese)
RadDennouAGCMFlux : | 放射フラックスの計算 |
RadDennouAGCMFinalize : | 終了処理 (モジュール内部の変数の割り付け解除) |
———— : | ———— |
RadDennouAGCMFlux : | Calculate radiation flux |
RadDennouAGCMFinalize : | Termination (deallocate variables in this module) |
Subroutine : |
リスタートファイルのクローズと, モジュール内部の変数の割り付け解除を行います.
Close a restart file, and deallocate variables in this module.
subroutine RadDennouAGCMFinalize ! ! リスタートファイルのクローズと, ! モジュール内部の変数の割り付け解除を行います. ! ! Close a restart file, and ! deallocate variables in this module. ! ! モジュール引用 ; USE statements ! ! リスタートデータ入出力 ! Restart data input/output ! use gtool_history, only: HistoryClose ! 宣言文 ; Declaration statements ! implicit none ! 実行文 ; Executable statement ! if ( .not. rad_DennouAGCM_inited ) return ! デフォルト値へ戻す ! Return to default values ! Old_Flux_saved = .false. ! リスタートファイルのクローズ ! close a restart file ! call HistoryClose( history = gthst_rst ) ! (inout) ! 割り付け解除 ! Deallocation ! if ( allocated( xy_IncomRadSFlux ) ) deallocate( xy_IncomRadSFlux ) if ( allocated( xy_InAngle ) ) deallocate( xy_InAngle ) if ( allocated( xy_TempSave ) ) deallocate( xy_TempSave ) if ( allocated( xyr_RadLFluxSave ) ) deallocate( xyr_RadLFluxSave ) if ( allocated( xyr_RadSFluxSave ) ) deallocate( xyr_RadSFluxSave ) if ( allocated( xyra_DelRadLFluxSave ) ) deallocate( xyra_DelRadLFluxSave ) rad_DennouAGCM_inited = .false. end subroutine RadDennouAGCMFinalize
Subroutine : | |||
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in)
| ||
xyz_QVap(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_SurfTemp(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xy_SurfAlbedo(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
| ||
flag_rst : | logical, intent(in), optional
|
温度, 比湿, 気圧から, 放射フラックスを計算します.
Calculate radiation flux from temperature, specific humidity, and air pressure.
subroutine RadDennouAGCMFlux( xyz_Temp, xyz_QVap, xyr_Press, xy_SurfTemp, xy_SurfAlbedo, xyr_RadLFlux, xyr_RadSFlux, xyra_DelRadLFlux, flag_rst ) ! ! 温度, 比湿, 気圧から, 放射フラックスを計算します. ! ! Calculate radiation flux from temperature, specific humidity, and ! air pressure. ! ! モジュール引用 ; USE statements ! ! 短波入射 (太陽入射) ! Short wave (insolation) incoming ! use rad_short_income, only: ShortIncoming ! 物理定数設定 ! Physical constants settings ! use constants, only: Grav ! $ g $ [m s-2]. ! 重力加速度. ! Gravitational acceleration ! 時刻管理 ! Time control ! use timeset, only: TimeN, EndTime, TimesetClockStart, TimesetClockStop ! リスタートデータ出力 ! Restart data output ! use gtool_history, only: HistoryPut, HistorySetTime ! デバッグ用ユーティリティ ! Utilities for debug ! use dc_trace, only: DbgMessage, BeginSub, EndSub ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyz_Temp (0:imax-1, 1:jmax, 1:kmax) ! $ T $ . 温度. Temperature real(DP), intent(in):: xyz_QVap (0:imax-1, 1:jmax, 1:kmax) ! $ q $ . 比湿. Specific humidity real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax) ! $ \hat{p} $ . 気圧 (半整数レベル). ! Air pressure (half level) real(DP), intent(in):: xy_SurfTemp (0:imax-1, 1:jmax) ! 地表面温度. ! Surface temperature real(DP), intent(in):: xy_SurfAlbedo (0:imax-1, 1:jmax) ! 地表アルベド. ! Surface albedo real(DP), intent(out):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Longwave flux real(DP), intent(out):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax) ! 短波 (日射) フラックス. ! Shortwave (insolation) flux real(DP), intent(out):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Long wave flux tendency with surface temperature (for save) logical, intent(in), optional:: flag_rst ! リスタートであることを示すフラグ. ! .true. が与えられる場合, ! 長波放射, 短波放射に関するリスタート ! ファイルが必要になります. ! リスタートファイルに関する情報は ! NAMELIST#rad_DennouAGCM_nml ! で指定されます. ! デフォルトは .false. です. ! ! Flag for restart. ! If .true. is given, ! a restart file for long radiation ! and short radiation is needed. ! Information about the restart file ! is specified by "NAMELIST#rad_DennouAGCM_nml". ! Default value is .false. ! ! 作業変数 ! Work variables ! real(DP):: xyr_ColDenQVap (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho q \, dz $ . ! 鉛直層 k より上空の水蒸気のカラム密度. ! Column density of water vapor above vertical level k. real(DP):: xyr_ColDenDryAir (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho \, dz $ . ! 鉛直層 k より上空の空気のカラム密度. ! Column density of air above vertical level k. integer:: k ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction logical:: flag_rst_output ! リスタートファイル出力のフラグ. ! Flag for output of a restart file real(DP) :: xy_CosSZA(0:imax-1, 1:jmax) real(DP) :: DistFromStarScld real(DP) :: DiurnalMeanFactor integer:: i integer:: j ! 実行文 ; Executable statement ! ! 計算時間計測開始 ! Start measurement of computation time ! call TimesetClockStart( module_name ) ! 初期化 ! Initialization ! if ( .not. rad_DennouAGCM_inited ) call RadInit( flag_rst ) ! 鉛直層 k より上空のカラム密度の計算 ! Calculate column density above vertical level k ! xyr_ColDenQVap (:,:,kmax) = 0. xyr_ColDenDryAir(:,:,kmax) = 0. do k = kmax-1, 0, -1 xyr_ColDenQVap(:,:,k) = xyr_ColDenQVap(:,:,k+1) + xyz_QVap(:,:,k+1) * ( xyr_Press(:,:,k) - xyr_Press(:,:,k+1) ) / Grav xyr_ColDenDryAir(:,:,k) = xyr_ColDenDryAir(:,:,k+1) + ( xyr_Press(:,:,k) - xyr_Press(:,:,k+1) ) / Grav end do ! 長波フラックスの算出 ! Calculate long wave flux ! if ( ( TimeN - PrevTimeLong >= IntTimeLong ) .or. ( .not. Old_Flux_saved ) ) then if ( .not. Old_Flux_saved ) then PrevTimeLong = TimeN else PrevTimeLong = PrevTimeLong + IntTimeLong end if call LongFlux( xyz_Temp, xy_SurfTemp, xyr_ColDenQVap, xyr_ColDenDryAir, xyr_RadLFlux, xyra_DelRadLFlux ) ! (out) ! 前回の値を利用 ! Use values in last time ! else xyr_RadLFlux = xyr_RadLFluxSave xyra_DelRadLFlux = xyra_DelRadLFluxSave if ( .not. flag_rst_input ) then do k = 0, kmax xyr_RadLFlux(:,:,k) = xyr_RadLFlux(:,:,k) + xyra_DelRadLFlux(:,:,k,1) * ( xyz_Temp(:,:,1) - xy_TempSave ) xyra_DelRadLFlux(:,:,k,1) = xyra_DelRadLFlux(:,:,k,1) / ( xy_TempSave**3 ) * ( xyz_Temp(:,:,1)**3 ) end do else flag_rst_input = .false. call DbgMessage( '%c: restart data is used. ', c1 = module_name ) end if end if ! 短波 (日射) フラックスの算出 ! Calculate short wave (insolation) ! if ( ( TimeN - PrevTimeShort >= IntTimeShort ) .or. ( .not. Old_Flux_saved ) ) then if ( .not. Old_Flux_saved ) then PrevTimeShort = TimeN else PrevTimeShort = PrevTimeShort + IntTimeShort end if ! 短波入射の計算 ! Calculate short wave (insolation) incoming radiation ! call ShortIncoming( xy_CosZet = xy_CosSZA, xy_InAngle = xy_InAngle, DistFromStarScld = DistFromStarScld, DiurnalMeanFactor = DiurnalMeanFactor ) do j = 1, jmax do i = 0, imax-1 if ( xy_CosSZA(i,j) > 0.0_DP ) then xy_IncomRadSFlux(i,j) = - ( SolarConst / DistFromStarScld** 2 ) * xy_CosSZA(i,j) else xy_IncomRadSFlux(i,j) = 0. end if end do end do ! Correction of incoming solar flux when diurnal mean insolation is assumed. ! If diurnal mean insolation is assumed, DiurnalMeanFactor = 1 / PI. ! If diurnal mean insolation is not assumed, DiurnalMeanFactor = 1. xy_IncomRadSFlux = xy_IncomRadSFlux * DiurnalMeanFactor ! 大気アルベドの考慮 ! Taking atmospheric albedo into consideration ! xy_IncomRadSFlux = xy_IncomRadSFlux * ( 1.0d0 - ShortAtmosAlbedo ) ! 短波フラックスの計算 ! Calculate short wave (insolation) flux ! xyr_RadSFlux(:,:,0:kmax-1) = 0. xyr_RadSFlux(:,:, kmax ) = xy_IncomRadSFlux call ShortFlux( xyr_ColDenQVap, xyr_ColDenDryAir, xy_SurfAlbedo, xyr_RadSFlux ) ! (inout) ! 前回の値を利用 ! Use values in last time ! else xyr_RadSFlux = xyr_RadSFluxSave end if ! 今回計算した値を保存 ! Save calculated values in this time ! xy_TempSave = xyz_Temp (:,:,1) xyr_RadLFluxSave = xyr_RadLFlux xyr_RadSFluxSave = xyr_RadSFlux xyra_DelRadLFluxSave = xyra_DelRadLFlux if ( .not. Old_Flux_saved ) Old_Flux_saved = .true. ! リスタートファイルの出力タイミングのチェック ! Check output timing of a restart file ! !! !! old code to be deleted !! ! flag_rst_output = ( TimeN - PrevRstOutputTime >= RstFileIntTime ) ! if ( TimeN >= EndTime .and. .not. flag_rst_output_end ) then ! flag_rst_output = .true. ! flag_rst_output_end = .true. ! end if ! flag_rst_output = ( .not. TimeN == PrevRstOutputTime ) .and. flag_rst_output if ( TimeN - PrevRstOutputTime >= RstFileIntTime ) then flag_rst_output = .true. else flag_rst_output = .false. end if if ( TimeN >= EndTime .and. .not. flag_rst_output_end ) then flag_rst_output = .true. flag_rst_output_end = .true. end if if ( ( .not. TimeN == PrevRstOutputTime ) .and. flag_rst_output ) then flag_rst_output = .true. else flag_rst_output = .false. end if if ( flag_rst_output ) then ! 次回用に, 今回の出力 (希望) 時刻 を保存 ! Save output time (expected) in this time, for next time ! PrevRstOutputTime = PrevRstOutputTime + RstFileIntTime ! 時刻の設定 ! Set time ! call HistorySetTime( timed = TimeN, history = gthst_rst ) ! データ出力 ! Data output ! call HistoryPut( 'PrevTimeLong', PrevTimeLong, history = gthst_rst ) ! (in) call HistoryPut( 'PrevTimeShort', PrevTimeShort, history = gthst_rst ) ! (in) call HistoryPut( 'SurfTemp', xy_TempSave, history = gthst_rst ) ! (in) call HistoryPut( 'RadLFlux', xyr_RadLFluxSave, history = gthst_rst ) ! (in) call HistoryPut( 'RadSFlux', xyr_RadSFluxSave, history = gthst_rst ) ! (in) call HistoryPut( 'DelRadLFlux', xyra_DelRadLFluxSave, history = gthst_rst ) ! (in) end if ! 計算時間計測一時停止 ! Pause measurement of computation time ! call TimesetClockStop( module_name ) end subroutine RadDennouAGCMFlux
Variable : | |||
rad_DennouAGCM_inited = .false. : | logical, save, public
|
Variable : | |||
DelTimeLongUnit : | character(STRING), save
|
Variable : | |||
DelTimeLongValue : | real(DP), save
|
Variable : | |||
DelTimeShortUnit : | character(STRING), save
|
Variable : | |||
DelTimeShortValue : | real(DP), save
|
Subroutine : |
依存モジュールの初期化チェック
Check initialization of dependency modules
subroutine InitCheck ! ! 依存モジュールの初期化チェック ! ! Check initialization of dependency modules ! モジュール引用 ; USE statements ! ! NAMELIST ファイル入力に関するユーティリティ ! Utilities for NAMELIST file input ! use namelist_util, only: namelist_util_inited ! 格子点設定 ! Grid points settings ! use gridset, only: gridset_inited ! 物理定数設定 ! Physical constants settings ! use constants, only: constants_inited ! 座標データ設定 ! Axes data settings ! use axesset, only: axesset_inited ! 時刻管理 ! Time control ! use timeset, only: timeset_inited ! リスタートデータ入出力 ! Restart data input/output ! use restart_file_io, only: restart_file_io_inited ! 実行文 ; Executable statement ! if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' ) if ( .not. gridset_inited ) call MessageNotify( 'E', module_name, '"gridset" module is not initialized.' ) if ( .not. constants_inited ) call MessageNotify( 'E', module_name, '"constants" module is not initialized.' ) if ( .not. axesset_inited ) call MessageNotify( 'E', module_name, '"axesset" module is not initialized.' ) if ( .not. timeset_inited ) call MessageNotify( 'E', module_name, '"timeset" module is not initialized.' ) if ( .not. restart_file_io_inited ) call MessageNotify( 'E', module_name, '"restart_file_io" module is not initialized.' ) end subroutine InitCheck
Variable : | |||
IntTimeLong : | real(DP), save
|
Variable : | |||
IntTimeShort : | real(DP), save
|
Variable : | |||
LongAbsorpCoefDryAir(1:MaxNmlArySize) : | real(DP), save
|
Variable : | |||
LongAbsorpCoefQVap(1:MaxNmlArySize) : | real(DP), save
|
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)
| ||
xyr_ColDenQVap(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_ColDenDryAir(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_RadLFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(out)
| ||
xyra_DelRadLFlux(0:imax-1, 1:jmax, 0:kmax, 0:1) : | real(DP), intent(out)
|
長波フラックスの計算
Calculate long wave flux
subroutine LongFlux( xyz_Temp, xy_SurfTemp, xyr_ColDenQVap, xyr_ColDenDryAir, xyr_RadLFlux, xyra_DelRadLFlux ) ! ! 長波フラックスの計算 ! ! Calculate long wave flux ! ! モジュール引用 ; USE statements ! ! 物理定数設定 ! Physical constants settings ! use constants, only: StB ! $ \sigma_{SB} $ . ! ステファンボルツマン定数. ! Stefan-Boltzmann constant ! 宣言文 ; Declaration statements ! implicit none 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):: xyr_ColDenQVap (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho q \, dz $ . ! 鉛直層 k より上空の水蒸気のカラム密度. ! Column density of water vapor above vertical level k. real(DP), intent(in):: xyr_ColDenDryAir (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho \, dz $ . ! 鉛直層 k より上空の空気のカラム密度. ! Column density of air above vertical level k. real(DP), intent(out):: xyr_RadLFlux (0:imax-1, 1:jmax, 0:kmax) ! 長波フラックス. ! Longwave flux real(DP), intent(out):: xyra_DelRadLFlux (0:imax-1, 1:jmax, 0:kmax, 0:1) ! 長波地表温度変化. ! Surface temperature tendency with longwave ! 作業変数 ! Work variables ! real(DP):: xyrr_Trans(0:imax-1, 1:jmax, 0:kmax, 0:kmax) ! 透過係数. ! Transmission coefficient real(DP):: xyz_PiB (0:imax-1, 1:jmax, 1:kmax) ! $ \pi B = \sigma T^{4} $ real(DP):: xy_SurfPiB (0:imax-1, 1:jmax) ! 地表面の $ \pi B $ . ! $ \pi B $ on surface real(DP):: xy_PiDBDT1 (0:imax-1, 1:jmax) ! $ \pi DBDT = 4 \sigma T^{3} at lowest level$ real(DP):: xy_SurfPiDBDT (0:imax-1, 1:jmax) ! 地表面の $ \pi DBDT $ . ! $ \pi B $ on surface real(DP):: BandWeightSum ! バンドウェイトの和 ! Sum of band weights integer:: k, kk ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction integer:: bn ! 波長について回る DO ループ用作業変数 ! Work variables for DO loop in wavenumber bands ! 実行文 ; Executable statement ! ! バンドウェイトの設定 ! Configure band weight ! BandWeightSum = 0. do bn = 1, LongBandNum BandWeightSum = BandWeightSum + LongBandWeight(bn) end do do bn = 1, LongBandNum LongBandWeight(bn) = LongBandWeight(bn) / BandWeightSum end do ! 透過関数計算 ! Calculate transmission functions ! ! Initialization ! xyrr_Trans = 0. ! do bn = 1, LongBandNum do k = 0, kmax do kk = k, kmax xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,k,kk) + LongBandWeight(bn) * exp( - LongPathLengthFact * ( LongAbsorpCoefQVap(bn) * abs( xyr_ColDenQVap(:,:,kk) - xyr_ColDenQVap(:,:,k) ) + LongAbsorpCoefDryAir(bn) * abs( xyr_ColDenDryAir(:,:,kk) - xyr_ColDenDryAir(:,:,k) ) ) ) end do end do end do do k = 0, kmax do kk = 0, k-1 xyrr_Trans(:,:,k,kk) = xyrr_Trans(:,:,kk,k) end do end do ! $ \pi B $, $ \pi DBDT $ の計算 ! Calculate $ \pi B $ and $ \pi DBDT $ ! xyz_PiB = StB * ( xyz_Temp**4 ) xy_SurfPiB = StB * ( xy_SurfTemp**4 ) xy_PiDBDT1 = 4.0_DP * xyz_PiB(:,:,1) / xyz_Temp(:,:,1) xy_SurfPiDBDT = 4.0_DP * xy_SurfPiB / xy_SurfTemp ! Integrate radiative transfer equation ! do k = 0, kmax xyr_RadLFlux(:,:,k) = xy_SurfPiB * xyrr_Trans(:,:,k,0) do kk = 0, kmax-1 xyr_RadLFlux(:,:,k) = xyr_RadLFlux(:,:,k) - xyz_PiB(:,:,kk+1) * ( xyrr_Trans(:,:,k,kk) - xyrr_Trans(:,:,k,kk+1) ) end do end do ! 長波地表温度変化の計算 ! Calclate surface temperature tendency with long wave ! do k = 0, kmax xyra_DelRadLFlux(:,:,k,0) = xy_SurfPiDBDT * xyrr_Trans(:,:,k,0) xyra_DelRadLFlux(:,:,k,1) = xy_PiDBDT1 * ( xyrr_Trans(:,:,k,1) - xyrr_Trans(:,:,k,0) ) end do end subroutine LongFlux
Variable : | |||
LongPathLengthFact : | real(DP), save
|
Variable : | |||
Old_Flux_saved = .false. : | logical, save
|
Variable : | |||
PrevRstOutputTime : | real(DP), save
|
Variable : | |||
PrevTimeLong : | real(DP), save
|
Variable : | |||
PrevTimeShort : | real(DP), save
|
Subroutine : | |||
flag_rst : | logical, intent(in), optional
|
rad_DennouAGCM モジュールの初期化を行います. NAMELIST#rad_DennouAGCM_nml の読み込みはこの手続きで行われます.
"rad_DennouAGCM" module is initialized. "NAMELIST#rad_DennouAGCM_nml" is loaded in this procedure.
This procedure input/output NAMELIST#rad_DennouAGCM_nml .
subroutine RadInit( flag_rst ) ! ! rad_DennouAGCM モジュールの初期化を行います. ! NAMELIST#rad_DennouAGCM_nml の読み込みはこの手続きで行われます. ! ! "rad_DennouAGCM" module is initialized. ! "NAMELIST#rad_DennouAGCM_nml" is loaded in this procedure. ! ! モジュール引用 ; USE statements ! ! 出力ファイルの基本情報 ! Basic information for output files ! use fileset, only: FileTitle, FileSource, FileInstitution ! データファイルを最終的に変更した組織/個人. ! Institution or person that changes data files for the last time ! 物理定数設定 ! Physical constants settings ! use constants, only: PI ! $ \pi $ . ! 円周率. Circular constant ! 座標データ設定 ! Axes data settings ! use axesset, only: x_Lon, x_Lon_Weight, y_Lat, y_Lat_Weight, z_Sigma, r_Sigma, z_DelSigma ! $ \Delta \sigma $ (整数). ! $ \Delta \sigma $ (Full) ! 時刻管理 ! Time control ! use timeset, only: RestartTime ! リスタート開始時刻. ! Retart time of calculation ! NAMELIST ファイル入力に関するユーティリティ ! Utilities for NAMELIST file input ! use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid ! 暦と日時の取り扱い ! Calendar and Date handler ! use dc_calendar, only: DCCalConvertByUnit ! ファイル入出力補助 ! File I/O support ! use dc_iounit, only: FileOpen ! 種別型パラメタ ! Kind type parameter ! use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output ! 組み込み関数 PRESENT の拡張版関数 ! Extended functions of intrinsic function "PRESENT" ! use dc_present, only: present_and_true ! リスタートデータ入出力 ! Restart data input/output ! use gtool_history, only: HistoryCreate, HistoryAddAttr, HistoryAddVariable, HistoryPut, HistoryGet, HistoryGetAttr ! 文字列操作 ! Character handling ! use dc_string, only: toChar ! 宣言文 ; Declaration statements ! implicit none logical, intent(in), optional:: flag_rst ! リスタートであることを示すフラグ. ! .true. が与えられる場合, ! 長波放射, 短波放射に関するリスタート ! ファイルが必要になります. ! リスタートファイルに関する情報は ! NAMELIST#rad_DennouAGCM_nml ! で指定されます. ! デフォルトは .false. です. ! ! Flag for restart. ! If .true. is given, ! a restart file for long radiation ! and short radiation is needed. ! Information about the restart file ! is specified by "NAMELIST#rad_DennouAGCM_nml". ! Default value is .false. ! character(STRING):: RstInputFile ! 入力するリスタートデータのファイル名 ! Filename of input restart data character(STRING):: RstOutputFile ! 出力するリスタートデータのファイル名 ! Filename of output restart data character(STRING):: time_range ! 時刻の指定. ! Specification of time character(TOKEN):: dummy_str ! 入力チェック用のダミー変数 ! Dummy variable for check of input logical:: get_err ! 入力時のエラーフラグ. ! Error flag for input integer:: unit_nml ! NAMELIST ファイルオープン用装置番号. ! Unit number for NAMELIST file open integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT. ! IOSTAT of NAMELIST read character(STRING):: title_msg ! 表題に付加するメッセージ. ! Message added to title real(DP):: origin_time ! 計算開始時刻. ! Start time of calculation character(12):: time_unit ! 日時の単位. Units of date and time logical:: flag_mpi_init ! NAMELIST 変数群 ! NAMELIST group name ! namelist /rad_DennouAGCM_nml/ SolarConst, DelTimeLongValue, DelTimeLongUnit, DelTimeShortValue, DelTimeShortUnit, LongBandNum, LongAbsorpCoefQVap, LongAbsorpCoefDryAir, LongBandWeight, LongPathLengthFact, ShortBandNum, ShortAbsorpCoefQVap, ShortAbsorpCoefDryAir, ShortBandWeight, ShortSecScat, ShortAtmosAlbedo, RstInputFile, RstOutputFile ! ! デフォルト値については初期化手続 "rad_DennouAGCM#RadInit" ! のソースコードを参照のこと. ! ! Refer to source codes in the initialization procedure ! "rad_DennouAGCM#RadInit" for the default values. ! ! 実行文 ; Executable statement ! if ( rad_DennouAGCM_inited ) return call InitCheck ! flag_mpi_init = .false. flag_mpi_init = .true. ! デフォルト値の設定 ! Default values settings ! ! 長波フラックス用情報 ! Information for long wave flux ! PrevTimeLong = RestartTime DelTimeLongValue = 3.0_DP DelTimeLongUnit = 'hrs.' LongBandNum = 4 LongAbsorpCoefQVap = -999.9_DP LongAbsorpCoefDryAir = -999.9_DP LongBandWeight = -999.9_DP LongAbsorpCoefQVap (1:LongBandNum) = (/ 8.0_DP, 1.0_DP, 0.1_DP, 0.0_DP /) LongAbsorpCoefDryAir (1:LongBandNum) = (/ 0.0_DP, 0.0_DP, 0.0_DP, 5.0e-5_DP /) LongBandWeight (1:LongBandNum) = (/ 0.2_DP, 0.1_DP, 0.1_DP, 0.6_DP /) LongPathLengthFact = 1.5_DP ! 短波フラックス用情報 ! Information for short wave flux ! SolarConst = 1380.0_DP PrevTimeShort = RestartTime DelTimeShortValue = 1.0_DP DelTimeShortUnit = 'hrs.' ShortBandNum = 1 ShortAbsorpCoefQVap = -999.9_DP ShortAbsorpCoefDryAir = -999.9_DP ShortBandWeight = -999.9_DP ShortAbsorpCoefQVap (1:ShortBandNum) = (/ 0.002_DP /) ShortAbsorpCoefDryAir (1:ShortBandNum) = (/ 0.0_DP /) ShortBandWeight (1:ShortBandNum) = (/ 1.0_DP /) ShortSecScat = 1.66_DP ShortAtmosAlbedo = 0.2_DP ! リスタートファイル情報 ! Information about a restart file ! RstInputFile = '' RstOutputFile = 'rst_rad.nc' ! 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_DennouAGCM_nml, iostat = iostat_nml ) ! (out) close( unit_nml ) call NmlutilMsg( iostat_nml, module_name ) ! (in) end if ! 時間間隔の処理 ! Handle interval time ! IntTimeLong = DCCalConvertByUnit( DelTimeLongValue, DelTimeLongUnit, 'sec' ) ! (in) IntTimeShort = DCCalConvertByUnit( DelTimeShortValue, DelTimeShortUnit, 'sec' ) ! (in) ! バンド数, 吸収係数, バンドウェイトのチェック ! Check number of band, absorption coefficients, band weight ! call NmlutilAryValid( module_name, LongAbsorpCoefQVap, 'LongAbsorpCoefQVap', LongBandNum, 'LongBandNum' ) ! (in) call NmlutilAryValid( module_name, LongAbsorpCoefDryAir, 'LongAbsorpCoefDryAir', LongBandNum, 'LongBandNum' ) ! (in) call NmlutilAryValid( module_name, LongBandWeight, 'LongBandWeight', LongBandNum, 'LongBandNum' ) ! (in) call NmlutilAryValid( module_name, ShortAbsorpCoefQVap, 'ShortAbsorpCoefQVap', ShortBandNum, 'ShortBandNum' ) ! (in) call NmlutilAryValid( module_name, ShortAbsorpCoefDryAir, 'ShortAbsorpCoefDryAir', ShortBandNum, 'ShortBandNum' ) ! (in) call NmlutilAryValid( module_name, ShortBandWeight, 'ShortBandWeight', ShortBandNum, 'ShortBandNum' ) ! (in) ! 短波入射用変数の割付 ! Allocate variables for short wave (insolation) incoming radiation ! allocate( xy_IncomRadSFlux (0:imax-1, 1:jmax) ) allocate( xy_InAngle (0:imax-1, 1:jmax) ) ! 保存用の変数の割り付け ! Allocate variables for saving ! allocate( xy_TempSave (0:imax-1, 1:jmax) ) allocate( xyr_RadLFluxSave (0:imax-1, 1:jmax, 0:kmax) ) allocate( xyr_RadSFluxSave (0:imax-1, 1:jmax, 0:kmax) ) allocate( xyra_DelRadLFluxSave (0:imax-1, 1:jmax, 0:kmax, 0:1) ) ! リスタートファイルの入力 ! Input restart file ! if ( present_and_true( flag_rst ) ) then if ( trim(RstInputFile) == '' ) then call MessageNotify( 'E', module_name, 'a restart file is needed. ' // 'Specify the restart file to "RstInputFile" in NAMELIST "rad_DennouAGCM_nml"' ) end if ! 時刻情報の取得 ! Get time information ! time_range = 'time=' // toChar( RestartTime ) ! ファイルの有無を確認 ! Conform an existence of an input file ! call HistoryGetAttr( RstInputFile, 'lon', 'units', dummy_str, flag_mpi_split = flag_mpi_init, err = get_err ) ! (out) if ( get_err ) then call MessageNotify( 'E', module_name, 'restart/initial data file "%c" is not found.', c1 = trim(RstInputFile) ) end if ! 入力 ! Input ! call HistoryGet( RstInputFile, 'SurfTemp', xy_TempSave, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional call HistoryGet( RstInputFile, 'RadLFlux', xyr_RadLFluxSave, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional call HistoryGet( RstInputFile, 'RadSFlux', xyr_RadSFluxSave, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional call HistoryGet( RstInputFile, 'DelRadLFlux', xyra_DelRadLFluxSave, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional call HistoryGet( RstInputFile, 'PrevTimeLong', PrevTimeLong, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional call HistoryGet( RstInputFile, 'PrevTimeShort', PrevTimeShort, range = time_range, flag_mpi_split = flag_mpi_init ) ! (in) optional flag_rst_input = .true. Old_Flux_saved = .true. else RstInputFile = '' flag_rst_input = .false. Old_Flux_saved = .false. end if ! 出力時間間隔の設定 ! Configure time interval of output ! PrevRstOutputTime = RestartTime ! リスタートファイルの作成 ! Create a restart file ! title_msg = ' restart data for "' // module_name // '" module' !!$ origin_time = EvalByUnit( RestartTime, RstFileIntUnit ) !!$ time_unit = RstFileIntUnit time_unit = 'sec' call HistoryCreate( file = RstOutputFile, title = trim(FileTitle) // trim(title_msg), source = FileSource, institution = FileInstitution, dims = (/ 'lon ', 'lat ', 'sig ', 'sigm ', 'sorbl', 'time ' /), dimsizes = (/ imax, jmax, kmax, kmax + 1, 2, 0 /), longnames = (/ 'longitude ', 'latitude ', 'sigma at layer midpoints ', 'sigma at layer end-points (half level)', 'surface or bottom layer ', 'time ' /), units = (/ 'degree_east ', 'degree_north', '1 ', '1 ', '1 ', time_unit /), xtypes = (/'double', 'double', 'double', 'double', 'int ', 'double'/), origind = RestartTime, intervald = RstFileIntValue, flag_mpi_split = flag_mpi_init, history = gthst_rst ) ! (out) optional ! 座標データの設定 ! Axes data settings ! call HistoryAddAttr( 'lon', attrname = 'standard_name', value = 'longitude', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'lat', attrname = 'standard_name', value = 'latitude', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'sig', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'sigm', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'time', attrname = 'standard_name', value = 'time', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'sig', attrname = 'positive', value = 'down', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'sigm', attrname = 'positive', value = 'down', history = gthst_rst ) ! (inout) call HistoryPut( 'lon', x_Lon / PI * 180.0_DP, history = gthst_rst ) ! (inout) call HistoryPut( 'lat', y_Lat / PI * 180.0_DP, history = gthst_rst ) ! (inout) call HistoryPut( 'sig', z_Sigma, history = gthst_rst ) ! (inout) call HistoryPut( 'sigm', r_Sigma, history = gthst_rst ) ! (inout) call HistoryPut( 'sorbl', (/ 0, 1 /), history = gthst_rst ) ! (inout) ! 座標重みの設定 ! Axes weights settings ! call HistoryAddVariable( 'lon_weight', (/'lon'/), 'weight for integration in longitude', 'radian', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'lon', attrname = 'gt_calc_weight', value = 'lon_weight', history = gthst_rst ) ! (inout) call HistoryPut( 'lon_weight', x_Lon_Weight, history = gthst_rst ) ! (inout) call HistoryAddVariable( 'lat_weight', (/'lat'/), 'weight for integration in latitude', units = 'radian', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'lat', attrname = 'gt_calc_weight', value = 'lat_weight', history = gthst_rst ) ! (inout) call HistoryPut( 'lat_weight', y_Lat_Weight, history = gthst_rst ) ! (inout) call HistoryAddVariable( 'sig_weight', (/'sig'/), 'weight for integration in sigma', '1', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddAttr( 'sig', attrname = 'gt_calc_weight', value = 'sig_weight', history = gthst_rst ) ! (inout) call HistoryPut( 'sig_weight', z_DelSigma, history = gthst_rst ) ! (inout) call HistoryAddVariable( 'PrevTimeLong', (/ 'time' /), 'previous time at which longwave flux is calculated', 'sec', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddVariable( 'PrevTimeShort', (/ 'time' /), 'previous time at which shortwave flux is calculated', 'sec', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddVariable( 'SurfTemp', (/ 'lon ', 'lat ', 'time' /), 'surface temperature', 'K', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddVariable( 'RadLFlux', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'longwave flux', 'W m-2', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddVariable( 'RadSFlux', (/ 'lon ', 'lat ', 'sigm', 'time' /), 'shortwave flux', 'W m-2', xtype = 'double', history = gthst_rst ) ! (inout) call HistoryAddVariable( 'DelRadLFlux', (/ 'lon ', 'lat ', 'sigm ', 'sorbl', 'time ' /), 'longwave flux tendency with surface temperature', 'W m-2 K-1', xtype = 'double', history = gthst_rst ) ! (inout) ! 印字 ; Print ! call MessageNotify( 'M', module_name, '----- Initialization Messages -----' ) call MessageNotify( 'M', module_name, 'Restart:' ) if ( trim(RstInputFile) == '' ) then call MessageNotify( 'M', module_name, ' InputFile = <no input>', c1 = trim( RstInputFile ) ) else call MessageNotify( 'M', module_name, ' InputFile = %c', c1 = trim( RstInputFile ) ) call MessageNotify( 'M', module_name, ' PrevTimeLong = %f [%c]', d = (/ PrevTimeLong /), c1 = 'sec' ) call MessageNotify( 'M', module_name, ' PrevTimeShort = %f [%c]', d = (/ PrevTimeShort /), c1 = 'sec' ) end if call MessageNotify( 'M', module_name, ' OutputFile = %c', c1 = trim( RstOutputFile ) ) call MessageNotify( 'M', module_name, ' IntTime = %f [%c] (same as IntTime in "restart_file_io" module)', d = (/ RstFileIntValue /), c1 = trim( RstFileIntUnit ) ) ! call MessageNotify( 'M', module_name, ' SolarConst = %f', d = (/ SolarConst /) ) call MessageNotify( 'M', module_name, 'DelTime:' ) call MessageNotify( 'M', module_name, ' DelTimeLong = %f [%c]', d = (/ DelTimeLongValue /), c1 = trim( DelTimeLongUnit ) ) call MessageNotify( 'M', module_name, ' DelTimeShort = %f [%c]', d = (/ DelTimeShortValue /), c1 = trim( DelTimeShortUnit ) ) ! call MessageNotify( 'M', module_name, 'LongFlux:' ) call MessageNotify( 'M', module_name, ' LongBandNum = %d', i = (/ LongBandNum /) ) call MessageNotify( 'M', module_name, ' LongAbsorpCoefQVap = (/ %*r /)', r = real( LongAbsorpCoefQVap(1:LongBandNum) ), n = (/ LongBandNum /) ) call MessageNotify( 'M', module_name, ' LongAbsorpCoefDryAir = (/ %*r /)', r = real( LongAbsorpCoefDryAir(1:LongBandNum) ), n = (/ LongBandNum /) ) call MessageNotify( 'M', module_name, ' LongBandWeight = (/ %*r /)', r = real( LongBandWeight(1:LongBandNum) ), n = (/ LongBandNum /) ) call MessageNotify( 'M', module_name, ' LongPathLengthFact = %f', d = (/ LongPathLengthFact /) ) ! call MessageNotify( 'M', module_name, 'ShortFlux:' ) call MessageNotify( 'M', module_name, ' ShortBandNum = %d', i = (/ ShortBandNum /) ) call MessageNotify( 'M', module_name, ' ShortAbsorpCoefQVap = (/ %*r /)', r = real( ShortAbsorpCoefQVap(1:ShortBandNum) ), n = (/ ShortBandNum /) ) call MessageNotify( 'M', module_name, ' ShortAbsorpCoefDryAir = (/ %*r /)', r = real( ShortAbsorpCoefDryAir(1:ShortBandNum) ), n = (/ ShortBandNum /) ) call MessageNotify( 'M', module_name, ' ShortBandWeight = (/ %*r /)', r = real( ShortBandWeight(1:ShortBandNum) ), n = (/ ShortBandNum /) ) call MessageNotify( 'M', module_name, ' ShortSecScat = %f', d = (/ ShortSecScat /) ) call MessageNotify( 'M', module_name, ' ShortAtmosAlbedo = %f', d = (/ ShortAtmosAlbedo /) ) call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) ) rad_DennouAGCM_inited = .true. end subroutine RadInit
Variable : | |||
ShortAbsorpCoefDryAir(1:MaxNmlArySize) : | real(DP), save
|
Variable : | |||
ShortAbsorpCoefQVap(1:MaxNmlArySize) : | real(DP), save
|
Subroutine : | |||
xyr_ColDenQVap(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xyr_ColDenDryAir(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(in)
| ||
xy_SurfAlbedo(0:imax-1, 1:jmax) : | real(DP), intent(in)
| ||
xyr_RadSFlux(0:imax-1, 1:jmax, 0:kmax) : | real(DP), intent(inout)
|
短波フラックスを計算します.
Calculate short wave flux.
subroutine ShortFlux( xyr_ColDenQVap, xyr_ColDenDryAir, xy_SurfAlbedo, xyr_RadSFlux ) ! ! 短波フラックスを計算します. ! ! Calculate short wave flux. ! ! モジュール引用 ; USE statements ! ! 宣言文 ; Declaration statements ! implicit none real(DP), intent(in):: xyr_ColDenQVap (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho q \, dz $ . ! 鉛直層 k より上空の水蒸気のカラム密度. ! Column density of water vapor above vertical level k. real(DP), intent(in):: xyr_ColDenDryAir (0:imax-1, 1:jmax, 0:kmax) ! $ \int_z^{\infty} \rho \, dz $ . ! 鉛直層 k より上空の空気のカラム密度. ! Column density of air above vertical level k. real(DP), intent(in):: xy_SurfAlbedo (0:imax-1, 1:jmax) ! 地表アルベド. ! Surface albedo real(DP), intent(inout):: xyr_RadSFlux (0:imax-1, 1:jmax, 0:kmax) ! 短波 (日射) フラックス. ! Shortwave (insolation) flux ! 作業変数 ! Work variables ! real(DP):: BandWeightSum ! バンドウェイトの和 ! Sum of band weights integer:: k ! 鉛直方向に回る DO ループ用作業変数 ! Work variables for DO loop in vertical direction integer:: bn ! 波長について回る DO ループ用作業変数 ! Work variables for DO loop in wavenumber bands ! 実行文 ; Executable statement ! ! バンドウェイトの設定 ! Configure band weight ! BandWeightSum = 0. do bn = 1, ShortBandNum BandWeightSum = BandWeightSum + ShortBandWeight(bn) end do do bn = 1, ShortBandNum ShortBandWeight(bn) = ShortBandWeight(bn) / BandWeightSum end do do bn = 1, ShortBandNum do k = 0, kmax ! 各レベルでの下向き透過 ! Downward transmission on each level ! if ( k /= kmax ) then xyr_RadSFlux(:,:,k) = xyr_RadSFlux(:,:,k) + ShortBandWeight(bn) * xyr_RadSFlux(:,:,kmax) * exp( - xy_InAngle(:,:) * ( ShortAbsorpCoefQVap(bn) * xyr_ColDenQVap(:,:,k) + ShortAbsorpCoefDryAir(bn) * xyr_ColDenDryAir(:,:,k) ) ) end if ! 各レベルでの上向き透過 ! Upward transmission on each level ! xyr_RadSFlux(:,:,k) = xyr_RadSFlux(:,:,k) - ShortBandWeight(bn) * xyr_RadSFlux(:,:,kmax) * exp( - xy_InAngle(:,:) * ( ShortAbsorpCoefQVap(bn) * xyr_ColDenQVap(:,:,0) + ShortAbsorpCoefDryAir(bn) * xyr_ColDenDryAir(:,:,0) ) ) * xy_SurfAlbedo * exp( - ShortSecScat * ( ShortAbsorpCoefQVap(bn) * ( xyr_ColDenQVap(:,:,0) - xyr_ColDenQVap(:,:,k) ) + ShortAbsorpCoefDryAir(bn) * ( xyr_ColDenDryAir(:,:,0) - xyr_ColDenDryAir(:,:,k) ) ) ) end do end do ! 吸収なしの場合 ! In the case of no absorption ! if ( ShortBandNum == 0 ) then do k = 0, kmax xyr_RadSFlux(:,:,k) = ( 1.0_DP - xy_SurfAlbedo ) * xyr_RadSFlux(:,:,kmax) end do end if end subroutine ShortFlux
Variable : | |||
ShortSecScat : | real(DP), save
|
Variable : | |||
flag_rst_input = .false. : | logical, save
|
Variable : | |||
flag_rst_output_end : | logical, save
|
Variable : | |||
gthst_rst : | type(GT_HISTORY), save
|
Constant : | |||
module_name = ‘rad_DennouAGCM‘ : | character(*), parameter
|
Constant : | |||
version = ’$Name: dcpam5-20110407 $’ // ’$Id: rad_DennouAGCM.f90,v 1.1 2011-03-27 02:20:12 yot Exp $’ : | character(*), parameter
|
Variable : | |||
xy_InAngle(:,:) : | real(DP), allocatable, save
|
Variable : | |||
xy_IncomRadSFlux(:,:) : | real(DP), allocatable, save
|
Variable : | |||
xy_TempSave(:,:) : | real(DP), allocatable, save
|
Variable : | |||
xyr_RadLFluxSave(:,:,:) : | real(DP), allocatable, save
|
Variable : | |||
xyr_RadSFluxSave(:,:,:) : | real(DP), allocatable, save
|
Variable : | |||
xyra_DelRadLFluxSave(:,:,:,:) : | real(DP), allocatable, save
|