!--------------------------------------------------------------------- ! Copyright (C) GFD Dennou Club, 2004. All rights reserved. !--------------------------------------------------------------------- != Subroutine GetDisturbVar ! ! * Developer: SUGIYAMA Ko-ichiro (sugiyama@gfd-dennou.org) ! * Version: $Id: getdisturbvar.f90,v 1.2 2006-06-17 10:21:44 sugiyama Exp $ ! * Tag Name: $Name: arare4-20080627 $ ! * Change History: ! !== Overview ! !擾乱場の初期値を netCDF ファイルから取得する. ! !== Error Handling ! !== Known Bugs ! !== Note ! !== Future Plans ! ! * 時間刻みのチェックだけでなく, 解像度のチェックおよび物理定数のチェックも必要だろう. ! subroutine GetDisturbVar(disturbfile, ReStartTime, & & pz_VelX_b, xr_VelZ_b, xz_Exner_b, & & xz_PotTemp_b, xz_Km_b, xz_Kh_b, & & xza_MixRtV_b, xza_MixRtC_b, xza_MixRtR_b, & & pz_VelX_n, xr_VelZ_n, xz_Exner_n, & & xz_PotTemp_n, xz_Km_n, xz_Kh_n, & & xza_MixRtV_n, xza_MixRtC_n, xza_MixRtR_n ) ! !擾乱場の初期値を netCDF ファイルから取得する. ! ! Dependency use gt4_history use fileset, only: InitFile use dc_trace, only: BeginSub, EndSub use dc_message ,only: MessageNotify use timeset, only: DelTimeLong use gridset, only: DimXMin, DimXMax, DimZMin, DimZMax, SpcNum !暗黙の型宣言禁止 implicit none ! Input character(*), intent(in) :: disturbfile !擾乱場のファイル ! Output real(8), intent(out) :: ReStartTime(2) real(8), intent(out) :: pz_VelX_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: pz_VelX_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xr_VelZ_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xr_VelZ_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Exner_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Exner_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_PotTemp_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_PotTemp_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Km_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Km_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Kh_n(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xz_Kh_b(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(out) :: xza_MixRtV_n(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) real(8), intent(out) :: xza_MixRtV_b(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) real(8), intent(out) :: xza_MixRtC_n(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) real(8), intent(out) :: xza_MixRtC_b(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) real(8), intent(out) :: xza_MixRtR_n(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) real(8), intent(out) :: xza_MixRtR_b(DimXMin:DimXMax, DimZMin:DimZMax,SpcNum) ! Work real :: DelTime character(30) :: name !変数名 call BeginSub("GetDisturbVar", & & fmt="%c", & & c1="Input initial disturbance variables.") write(*,*) InitFile ! Get a Value from netCDF File ! name = "t" ! call HistoryGet( InitFile, name, Time ) ReStartTime = (/ 0.0d0 , 2.0d0 /) !時間刻みのチェック DelTime = ReStartTime(2) - ReStartTime(1) ! if ( DelTime /= real(DelTimeLong, 4) ) then ! write(*,*) DelTime, DelTimeLong ! call MessageNotify("Error", "getdistrbvar", & ! & "DelTime for ReStartFile is not the same as DelTimeLong") ! stop ! end if ! Get a Value from netCDF File name = "VelX_n" call HistoryGet( InitFile, name, pz_VelX_n ) name = "VelX_b" call HistoryGet( InitFile, name, pz_VelX_b ) name = "VelZ_n" call HistoryGet( InitFile, name, xr_VelZ_n ) name = "VelZ_b" call HistoryGet( InitFile, name, xr_VelZ_b ) name = "Exner_n" call HistoryGet( InitFile, name, xz_Exner_n ) name = "Exner_b" call HistoryGet( InitFile, name, xz_Exner_b ) name = "PotTemp_n" call HistoryGet( InitFile, name, xz_PotTemp_n ) name = "PotTemp_b" call HistoryGet( InitFile, name, xz_PotTemp_b ) name = "Km_n" call HistoryGet( InitFile, name, xz_Km_n ) name = "Km_b" call HistoryGet( InitFile, name, xz_Km_b ) name = "Kh_n" call HistoryGet( InitFile, name, xz_Kh_n ) name = "Kh_b" call HistoryGet( InitFile, name, xz_Kh_b ) name = "MixRtV_n" call HistoryGet( InitFile, name, xza_MixRtV_n(:,:,1) ) name = "MixRtV_b" call HistoryGet( InitFile, name, xza_MixRtV_b(:,:,1) ) name = "MixRtC_n" call HistoryGet( InitFile, name, xza_MixRtC_n(:,:,1) ) name = "MixRtC_b" call HistoryGet( InitFile, name, xza_MixRtC_b(:,:,1) ) name = "MixRtR_n" call HistoryGet( InitFile, name, xza_MixRtR_n(:,:,1) ) name = "MixRtR_b" call HistoryGet( InitFile, name, xza_MixRtR_b(:,:,1) ) call EndSub("GetDisturbVar") end subroutine GetDisturbVar