Class ReStartFileIO
In: io/restartfileio.f90

リスタート用の場の情報を netCDF ファイルに出力するためのルーチン

Methods

Included Modules

gtool_history dc_message dc_string gridset timeset fileset basicset

Public Instance methods

Subroutine :

リスタートファイルのクローズ

[Source]

  subroutine ReStartFile_Close
    !
    !リスタートファイルのクローズ
    !
    
    !モジュール読み込み
    use gtool_history
    
    !暗黙の型宣言禁止
    implicit none

    !ファイルを閉じる
    call HistoryClose(rstat)
    
  end subroutine ReStartFile_Close
Subroutine :
ReStartTime(2) :real(8), intent(out)
xz_PotTempB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xz_ExnerB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
pz_VelXB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xr_VelZB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xza_MixRtB(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(out)
xz_KmB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xz_KhB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xz_PotTempN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xz_ExnerN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
pz_VelXN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xr_VelZN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xza_MixRtN(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(out)
xz_KmN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)
xz_KhN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(out)

リスタートファイルから情報取得

[Source]

  subroutine ReStartFile_Get( ReStartTime, xz_PotTempB, xz_ExnerB, pz_VelXB, xr_VelZB, xza_MixRtB, xz_KmB, xz_KhB, xz_PotTempN, xz_ExnerN, pz_VelXN, xr_VelZN, xza_MixRtN, xz_KmN, xz_KhN  )
    !
    !リスタートファイルから情報取得
    !

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(8), intent(out) :: ReStartTime(2)
    real(8), intent(out) :: pz_VelXN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xr_VelZN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_ExnerN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_PotTempN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_KmN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_KhN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xza_MixRtN(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8), intent(out) :: pz_VelXB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xr_VelZB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_ExnerB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_PotTempB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_KmB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_KhB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xza_MixRtB(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8)              :: DelTime
    real(8)              :: Var2D(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)              :: Var3D(DimXMin:DimXMax, DimZMin:DimZMax, 2)
    real(8)              :: Var3Ds(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8)              :: Var4D(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum, 2)
    real(8)              :: xz_ExnerBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場のエクスナー関数
    real(8)              :: xz_DensBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の密度
    real(8)              :: xz_PotTempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の温位
    real(8)              :: xz_VelSoundBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の音速
    real(8)              :: xz_PressBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の圧力
    real(8)              :: xz_TempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の温度
    real(8)              :: xza_MixRtBasicZ(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
                                               !基本場の混合比
    real(8)              :: xz_EffMolWtBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
                                               !基本場の分子量効果
    character(30)        :: name               !変数名
    character(10)        :: step(2)
    integer              :: t


    !-------------------------------------------------------------
    ! Get a Value from netCDF File
    !-------------------------------------------------------------    
    do t = 1, 2
      step(t) = 't=^' // adjustl(toChar(t))   
    end do

    do t = 1, 2
      name = "t"
      call HistoryGet( InitFile, name, ReStartTime(t), step(t) )
    end do

    do t = 1, 2
      name = "VelX"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    pz_VelXB = Var3D(:,:,1)
    pz_VelXN = Var3D(:,:,2)

    do t = 1, 2    
      name = "VelZ"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xr_VelZB = Var3D(:,:,1)
    xr_VelZN = Var3D(:,:,2)

    do t = 1, 2    
      name = "Exner"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_ExnerB = Var3D(:,:,1)
    xz_ExnerN = Var3D(:,:,2)

    do t = 1, 2    
      name = "PotTemp"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_PotTempB = Var3D(:,:,1)
    xz_PotTempN = Var3D(:,:,2)
    
    do t = 1, 2  
      name = "Km"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_KmB = Var3D(:,:,1)
    xz_KmN = Var3D(:,:,2)
    
    do t = 1, 2  
      name = "Kh"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_KhB = Var3D(:,:,1)
    xz_KhN = Var3D(:,:,2)
      
    do t = 1, 2
      name = "MixRt"
      call HistoryGet( InitFile, name, Var4D(:,:,:,t), step(t) )    
    end do
    xza_MixRtB = Var4D(:,:,:,1)
    xza_MixRtN = Var4D(:,:,:,2)
    
        
    !-------------------------------------------------------------
    ! 基本場の取得
    !-------------------------------------------------------------
    name = "DensBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_DensBasicZ = Var2D
    
    name = "ExnerBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_ExnerBasicZ = Var2D
    
    name = "PotTempBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_PotTempBasicZ = Var2D
    
    name = "VelSoundBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_VelSoundBasicZ = Var2D
    
    name = "TempBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_TempBasicZ = Var2D
    
    name = "PressBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_PressBasicZ = Var2D
    
    name = "MixRtBasicZ"
    call HistoryGet( InitFile, name, Var3Ds )
    xza_MixRtBasicZ = Var3Ds
    
    name = "EffMolWtBasicZ"
    call HistoryGet( InitFile, name, Var2D )
    xz_EffMolWtBasicZ = Var2D

    !----------------------------------------------------------
    ! 時間刻みのチェック
    !----------------------------------------------------------
    DelTime = ReStartTime(2) - ReStartTime(1)
    if ( DelTime /= DelTimeLong ) then 
      call MessageNotify( "W", "ReStartFile_Get", "DelTime in InitFile is not the same as DelTimeLong")
      call MessageNotify( "M", "ReStartFile_Get", "DelTime=%d", d=(/DelTime/) )
      call MessageNotify( "M", "ReStartFile_Get", "DelTime=Long%d", d=(/DelTimeLong/) )
    end if
    
    !----------------------------------------------------------
    ! BasicSet モジュールに値を設定
    !----------------------------------------------------------
    call BasicSetArray_Init( xz_PressBasicZ,    xz_ExnerBasicZ, xz_TempBasicZ, xz_PotTempBasicZ,  xz_DensBasicZ,  xz_VelSoundBasicZ, xza_MixRtBasicZ, xz_EffMolWtBasicZ )
    
  end subroutine ReStartFile_Get
Subroutine :

リスタートファイルの書き出し

[Source]

  subroutine ReStartFile_Open(  )
    !
    !リスタートファイルの書き出し
    !
   
    use basicset, only: xz_ExnerBasicZ, xz_DensBasicZ, xz_PotTempBasicZ, xz_VelSoundBasicZ, xz_PressBasicZ, xz_TempBasicZ, xza_MixRtBasicZ, xz_EffMolWtBasicZ   !基本場の分子量効果

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8)            :: SpcID(SpcNum)
    integer            :: N, M
    integer            :: s    

    SpcID = 0.0d0
    do s = 1, SpcNum
      SpcID(s) = real( s, 8 )
    end do
    
    N = size(s_X, 1)
    M = size(s_Z, 1)
    
    !-------------------------------------------------------------    
    ! ヒストリー作成
    !-------------------------------------------------------------  
    call HistoryCreate( file = ReStartFile, title = exptitle, source = expsrc, institution = expinst, dims=(/'x','z','s','t'/), dimsizes=(/N, M, SpcNum, 0/), longnames=(/'X-coordinate', 'Z-coordinate', 'Species Num ', 'Time        '/), units=(/'m','m','1','t'/), xtypes=(/'double', 'double', 'double', 'double'/), origin=0.0, interval=1.0, history=rstat )
    
    !-------------------------------------------------------------  
    ! 変数出力
    !-------------------------------------------------------------
    call HistoryPut('x', s_X, rstat)
    call HistoryPut('z', s_Z, rstat)
    call HistoryPut('s', SpcID, rstat)
    
    !無次元圧力の基本場
    call HistoryAddVariable( varname='ExnerBasicZ', dims=(/'x','z'/), longname='nondimensional pressure', units='1', xtype='double', history=rstat )
    
    !温位の基本場
    call HistoryAddVariable( varname='PotTempBasicZ', dims=(/'x','z'/), longname='potential temperature', units='K', xtype='double', history=rstat )
    
    !密度の基本場
    call HistoryAddVariable( varname='DensBasicZ', dims=(/'x','z'/), longname='density', units='Kg.m-3', xtype='double', history=rstat )    

    !音波速度の基本場
    call HistoryAddVariable( varname='VelSoundBasicZ', dims=(/'x','z'/), longname='sound velocity', units='m.s-2', xtype='double', history=rstat )
    
    !温度の基本場
    call HistoryAddVariable( varname='TempBasicZ', dims=(/'x','z'/), longname='Temperature of basic state', units='K', xtype='double', history=rstat ) 
    
    !圧力の基本場
    call HistoryAddVariable( varname='PressBasicZ', dims=(/'x','z'/), longname='Pressure of basic state', units='Pa', xtype='double', history=rstat ) 
    
    !水蒸気混合比の基本場
    call HistoryAddVariable( varname='MixRtBasicZ', dims=(/'x','z','s'/), longname='Mixing ratio of Condensible volatiles', units='kg.kg-1', xtype='double', history=rstat ) 
    
    !分子量効果
    call HistoryAddVariable( varname='EffMolWtBasicZ', dims=(/'x','z'/), longname='Effect of Mole Weight', units='1', xtype='double', history=rstat ) 
    
    !無次元圧力
    call HistoryAddVariable( varname='Exner', dims=(/'x','z','t'/), longname='nondimensional pressure', units='1', xtype='double', history=rstat )
    
    !温位の擾乱
    call HistoryAddVariable( varname='PotTemp', dims=(/'x','z','t'/), longname='virtual potential temperature', units='K', xtype='double', history=rstat )
    
    !速度
    call HistoryAddVariable( varname='VelX', dims=(/'x','z','t'/), longname='zonal velocity', units='m.s-1', xtype='double', history=rstat )
    
    !速度
    call HistoryAddVariable( varname='VelZ', dims=(/'x','z','t'/), longname='vertical velocity', units='m.s-1', xtype='double', history=rstat )
    
    !渦粘性係数
    call HistoryAddVariable( varname='Km', dims=(/'x','z','t'/), longname='Km', units='m2.s-1', xtype='double', history=rstat )
    
    !渦粘性係数
    call HistoryAddVariable( varname='Kh', dims=(/'x','z','t'/), longname='Kh', units='m2.s-1', xtype='double', history=rstat )
    
    !混合比
    call HistoryAddVariable( varname='MixRt', dims=(/'x','z','s','t'/), longname='Mixing Ratio', units='kg.kg-1"', xtype='double', history=rstat )
    
    !-------------------------------------------------------------
    ! 基本場のファイル出力
    !-------------------------------------------------------------
    call HistoryPut( 'DensBasicZ',     xz_DensBasicZ     , rstat)
    call HistoryPut( 'ExnerBasicZ',    xz_ExnerBasicZ    , rstat)
    call HistoryPut( 'PotTempBasicZ',  xz_PotTempBasicZ  , rstat)
    call HistoryPut( 'VelSoundBasicZ', xz_VelSoundBasicZ , rstat)
    call HistoryPut( 'TempBasicZ',     xz_TempBasicZ     , rstat)
    call HistoryPut( 'PressBasicZ',    xz_PressBasicZ    , rstat)
    call HistoryPut( 'MixRtBasicZ',    xza_MixRtBasicZ   , rstat)
    call HistoryPut( 'EffMolWtBasicZ', xz_EffMolWtBasicZ , rstat)
    
  end subroutine ReStartFile_Open
Subroutine :
Time :real(8), intent(in)
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(in)
xz_Km(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)

リスタートファイルに予報変数を書き出す

[Source]

  subroutine ReStartFile_OutPut( Time, xz_PotTemp, xz_Exner, pz_VelX, xr_VelZ, xza_MixRt, xz_Km , xz_Kh )
    !
    !リスタートファイルに予報変数を書き出す
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in)  :: Time
    real(8), intent(in)  :: pz_VelX(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xr_VelZ(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_Km(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_Kh(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    
  
    !------------------------------------------------------------------
    ! ファイル出力
    !------------------------------------------------------------------
    call HistoryPut( 't',       Time       , rstat)
    call HistoryPut( 'VelX',    pz_VelX    , rstat)
    call HistoryPut( 'VelZ',    xr_VelZ    , rstat)
    call HistoryPut( 'Exner',   xz_Exner   , rstat) 
    call HistoryPut( 'PotTemp', xz_PotTemp , rstat)
    call HistoryPut( 'Km',      xz_Km      , rstat)
    call HistoryPut( 'Kh',      xz_Kh      , rstat)
    call HistoryPut( 'MixRt',   xza_MixRt  , rstat)    
    
  end subroutine ReStartFile_OutPut

[Validate]