Class | ReStartFileIO |
In: |
../src/io/restartfileio.f90
|
リスタート用の場の情報を netCDF ファイルに出力するためのルーチン
Subroutine : |
リスタートファイルのクローズ
subroutine ReStartFileIO_Finalize ! !リスタートファイルのクローズ ! !モジュール読み込み use gtool_history !暗黙の型宣言禁止 implicit none !ファイルを閉じる call HistoryClose(rstat, quiet=.true.) end subroutine ReStartFileIO_Finalize
Subroutine : |
リスタートファイルから情報取得
subroutine ReStartFileio_BZ_Get() ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP) :: Var3D (imin:imax,jmin:jmax,kmin:kmax) real(DP) :: Var4D (imin:imax,jmin:jmax,kmin:kmax,ncmax) real(DP) :: xyz_ExnerBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場のエクスナー関数 real(DP) :: xyz_DensBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の密度 real(DP) :: xyz_PTempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温位 real(DP) :: xyz_VelSoundBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の音速 real(DP) :: xyz_PressBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の圧力 real(DP) :: xyz_TempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温度 real(DP) :: xyzf_QMixBZ (imin:imax,jmin:jmax,kmin:kmax, ncmax) !基本場の混合比 real(DP) :: xyz_EffMolWtBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の分子量効果 character(STRING) :: name !変数名 !------------------------------------------------------------- ! 基本場の取得 ! name = "DensBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_DensBZ = Var3D name = "ExnerBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_ExnerBZ = Var3D name = "PTempBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_PTempBZ = Var3D name = "VelSoundBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_VelSoundBZ = Var3D name = "TempBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_TempBZ = Var3D name = "PressBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_PressBZ = Var3D name = "QMixBZ" call HistoryGet( InputFile, name, Var4D, flag_mpi_split = FLAG_LIB_MPI ) xyzf_QMixBZ = Var4D name = "EffMolWtBZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyz_EffMolWtBZ = Var3D !---------------------------------------------------------- ! BasicSet モジュールに値を設定 !---------------------------------------------------------- call basicset_init( xyz_PressBZ, xyz_ExnerBZ, xyz_TempBZ, xyz_PTempBZ, xyz_DensBZ, xyz_VelSoundBZ, xyzf_QMixBZ, xyz_EffMolWtBZ ) end subroutine ReStartFileio_BZ_Get
Subroutine : |
リスタートファイルの書き出し
This procedure input/output NAMELIST#restartfileio_nml .
subroutine ReStartFileio_Init ! !リスタートファイルの書き出し ! !暗黙の型宣言禁止 implicit none !変数定義 real(4) :: SpcID(ncmax) integer :: N, L, M integer :: s integer :: unit !装置番号 !NAMELIST から情報を取得 NAMELIST /restartfileio_nml/ InputFile, OutputFile call FileOpen(unit, file=namelist_filename, mode='r') read(unit, NML=restartfileio_nml) close(unit) !確認 if (myrank == 0) then call MessageNotify( "M", "restartfileioIO_init", "InputFile = %c", c1=trim(InputFile) ) call MessageNotify( "M", "restartfileioIO_init", "OutputFile = %c", c1=trim(OutputFile) ) end if SpcID = 0.0d0 do s = 1, ncmax SpcID(s) = real( s, 4 ) end do N = size(x_X, 1) L = size(y_Y, 1) M = size(z_Z, 1) !------------------------------------------------------------- ! ヒストリー作成 !------------------------------------------------------------- call HistoryCreate( file = Outputfile, title = filetitle, source = filesource, institution = FileInstitution, dims=(/'x','y','z','s','t'/), dimsizes=(/N, L, M, ncmax, 0/), longnames=(/'X-coordinate', 'Y-coordinate', 'Z-coordinate', 'Species Num ', 'Time '/), units=(/'m ','m ','m ','1 ','sec'/), xtypes=(/'double', 'double', 'double', 'double', 'double'/), flag_mpi_split = FLAG_LIB_MPI, origin=0.0, interval=1.0, history=rstat, quiet=.true. ) !------------------------------------------------------------- ! 変数出力 !------------------------------------------------------------- call HistoryPut('x', x_X, rstat ) call HistoryPut('y', y_Y, rstat ) call HistoryPut('z', z_Z, rstat ) call HistoryPut('s', real(SpcID, 4), rstat ) !無次元圧力の基本場 call HistoryAddVariable( varname='ExnerBZ', dims=(/'x','y','z'/), longname='nondimensional pressure', units='1', xtype='double', history=rstat ) !温位の基本場 call HistoryAddVariable( varname='PTempBZ', dims=(/'x','y','z'/), longname='potential temperature', units='K', xtype='double', history=rstat ) !密度の基本場 call HistoryAddVariable( varname='DensBZ', dims=(/'x','y','z'/), longname='density', units='Kg.m-3', xtype='double', history=rstat ) !音波速度の基本場 call HistoryAddVariable( varname='VelSoundBZ', dims=(/'x','y','z'/), longname='sound velocity', units='m.s-2', xtype='double', history=rstat ) !温度の基本場 call HistoryAddVariable( varname='TempBZ', dims=(/'x','y','z'/), longname='Temperature of basic state', units='K', xtype='double', history=rstat ) !圧力の基本場 call HistoryAddVariable( varname='PressBZ', dims=(/'x','y','z'/), longname='Pressure of basic state', units='Pa', xtype='double', history=rstat ) !水蒸気混合比の基本場 call HistoryAddVariable( varname='QMixBZ', dims=(/'x','y','z','s'/), longname='Mixing ratio of Condensible volatiles', units='kg.kg-1', xtype='double', history=rstat ) !分子量効果 call HistoryAddVariable( varname='EffMolWtBZ', dims=(/'x','y','z'/), longname='Effect of Mole Weight', units='1', xtype='double', history=rstat ) !分子量効果 call HistoryAddVariable( varname='HumBZ', dims=(/'x','y','z','s'/), longname='Humidity', units='1', xtype='double', history=rstat ) !速度 call HistoryAddVariable( varname='VelX', dims=(/'x','y','z','t'/), longname='zonal velocity', units='m.s-1', xtype='double', history=rstat ) !速度 call HistoryAddVariable( varname='VelY', dims=(/'x','y','z','t'/), longname='meridional velocity', units='m.s-1', xtype='double', history=rstat ) !速度 call HistoryAddVariable( varname='VelZ', dims=(/'x','y','z','t'/), longname='vertical velocity', units='m.s-1', xtype='double', history=rstat ) !無次元圧力 call HistoryAddVariable( varname='Exner', dims=(/'x','y','z','t'/), longname='nondimensional pressure', units='1', xtype='double', history=rstat ) !温位の擾乱 call HistoryAddVariable( varname='PTemp', dims=(/'x','y','z','t'/), longname='virtual potential temperature', units='K', xtype='double', history=rstat ) !渦粘性係数 call HistoryAddVariable( varname='Km', dims=(/'x','y','z','t'/), longname='Km', units='m2.s-1', xtype='double', history=rstat ) !渦粘性係数 call HistoryAddVariable( varname='Kh', dims=(/'x','y','z','t'/), longname='Kh', units='m2.s-1', xtype='double', history=rstat ) !雲密度 call HistoryAddVariable( varname='CDens', dims=(/'x','y','z','t'/), longname='CDens', units='m2.s-1', xtype='double', history=rstat ) !混合比 call HistoryAddVariable( varname='QMix', dims=(/'x','y','z','s','t'/), longname='Mixing Ratio', units='kg.kg-1"', xtype='double', history=rstat ) end subroutine ReStartFileio_Init
Subroutine : | |
pyz_VelXB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
pyz_VelXN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xqz_VelYB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xqz_VelYN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyr_VelZB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyr_VelZN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_PTempB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_PTempN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_ExnerB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_ExnerN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyzf_QMixB(imin:imax,jmin:jmax,kmin:kmax,ncmax) : | real(DP), intent(out) |
xyzf_QMixN(imin:imax,jmin:jmax,kmin:kmax,ncmax) : | real(DP), intent(out) |
xyz_KmB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_KmN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_KhB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_KhN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_CDensB(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_CDensN(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
リスタートファイルから情報取得
subroutine ReStartFileio_Var_Get( pyz_VelXB, pyz_VelXN, xqz_VelYB, xqz_VelYN, xyr_VelZB, xyr_VelZN, xyz_PTempB, xyz_PTempN, xyz_ExnerB, xyz_ExnerN, xyzf_QMixB, xyzf_QMixN, xyz_KmB, xyz_KmN, xyz_KhB, xyz_KhN, xyz_CDensB, xyz_CDensN ) ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP), intent(out) :: pyz_VelXN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xqz_VelYN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyr_VelZN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_ExnerN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_PTempN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_KmN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_KhN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_CDensN (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyzf_QMixN (imin:imax,jmin:jmax,kmin:kmax,ncmax) real(DP), intent(out) :: pyz_VelXB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xqz_VelYB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyr_VelZB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_ExnerB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_PTempB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_KmB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_KhB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_CDensB (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyzf_QMixB (imin:imax,jmin:jmax,kmin:kmax,ncmax) character(STRING) :: name !変数名 character(STRING) :: TimeN = "" character(STRING) :: TimeB = "" !------------------------------------------------------------- ! Get a Value from netCDF File !------------------------------------------------------------- if (RestartTime /= 0.0d0) then timeB = 't=' // toChar( RestartTime - DelTimeLong) timeN = 't=' // toChar( RestartTime ) else timeB = 't=' // toChar( RestartTime ) timeN = 't=' // toChar( RestartTime ) end if ! write(*,*) TimeN, TimeB name = "PTemp" call HistoryGet( InputFile, name, xyz_PTempB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyz_PTempN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "VelX" call HistoryGet( InputFile, name, pyz_VelXB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, pyz_VelXN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "VelY" call HistoryGet( InputFile, name, xqz_VelYB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xqz_VelYN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "VelZ" call HistoryGet( InputFile, name, xyr_VelZB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyr_VelZN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "Exner" call HistoryGet( InputFile, name, xyz_ExnerB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyz_ExnerN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "Km" call HistoryGet( InputFile, name, xyz_KmB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyz_KmN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "Kh" call HistoryGet( InputFile, name, xyz_KhB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyz_KhN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "CDens" call HistoryGet( InputFile, name, xyz_CDensB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyz_CDensN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) name = "QMix" call HistoryGet( InputFile, name, xyzf_QMixB, range=TimeB, flag_mpi_split = FLAG_LIB_MPI ) call HistoryGet( InputFile, name, xyzf_QMixN, range=TimeN, flag_mpi_split = FLAG_LIB_MPI ) end subroutine ReStartFileio_Var_Get