!= Module ReStartFileIO
!
! Authors::   SUGIYAMA Ko-ichiro
! Version::   $Id: restartfileio.f90,v 1.13 2007/08/24 06:42:55 sugiyama Exp $ 
! Tag Name::  $Name: arare4-20071012 $
! Copyright:: Copyright (C) GFD Dennou Club, 2006. All rights reserved.
! License::   See COPYRIGHT[link:../../COPYRIGHT]
!
!== Overview 
!
!リスタート用の場の情報を netCDF ファイルに出力するためのルーチン
!
!== Error Handling
!
!== Known Bugs
!
!== Note
!
!速度を評価する格子位置と, 軸として入力してある座標値は整合的でないことに注意
!
!== Future Plans
!

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

  !モジュール読み込み
  use gt4_history, only: HistoryCreate, HistoryPut, HistoryGet, &
    &                    HistoryAddVariable, HistoryClose, gt_history
  use dc_message, only: MessageNotify
  use dc_string
  use gridset,  only: s_X,               &!X 座標軸(スカラー格子点)
    &                 s_Z,               &!Z 座標軸(スカラー格子点)
    &                 DimXMin,           &!X 方向の配列の下限
    &                 DimXMax,           &!X 方向の配列の上限
    &                 DimZMin,           &!Z 方向の配列の下限
    &                 DimZMax,           &!Z 方向の配列の上限
    &                 RegXMin,           &! 物理領域の X 方向の下限
    &                 RegXMax,           &! 物理領域の X 方向の上限
    &                 SpcNum              !凝縮成分の数
  use timeset,  only: DelTimeLong
  use fileset,  only: RestartFile,       &!リスタートファイル名
    &                 InitFile,          &!初期ファイル名
    &                 exptitle,          &!データの表題
    &                 expsrc,            &!データを作成する手順
    &                 expinst             !最終変更者・組織
  use basicset, only: BasicSetArray_Init, &!
    &                 PressBasis,    &!温位の基準圧力
    &                 GasRDry,       &!乾燥成分の定圧比熱
    &                 CpDry,         &!乾燥成分の定圧比熱
    &                 CvDry,         &!乾燥成分の定積比熱
    &                 MolWtDry,      &!乾燥成分の分子量
    &                 Grav,          &!重力加速度
    &                 SpcWetMolFr,   &!凝縮成分の初期モル比
    &                 MolWtWet,      &!凝縮成分の分子量
    &                 GasRUniv        !

  !暗黙の型宣言禁止
  implicit none
  
  !属性の指定
  private

  !関数を public に指定
  public ReStartFile_Open
  public ReStartFile_OutPut
  public ReStartFile_Close
  public ReStartFile_Get

  type(GT_HISTORY)   :: rstat

  save rstat 

contains

  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', 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'/),     &
      & varname='PotTempDist', 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',                                 &
      & xtype='double', history=rstat )
    
    !速度
    call HistoryAddVariable(                         &
      & varname='VelZ', dims=(/'x','z','t'/),        &
      & longname='vertical velocity',                &
      & units='m/s',                                 &
      & xtype='double', history=rstat )
    
    !渦粘性係数
    call HistoryAddVariable(                         &
      & varname='Km', dims=(/'x','z','t'/),          &
      & longname='Km',                               &
      & units='1',                                   &
      & xtype='double', history=rstat )
    
    !渦粘性係数
    call HistoryAddVariable(                         &
      & varname='Kh', dims=(/'x','z','t'/),          &
      & longname='Kh',                               &
      & units='1',                                   &
      & xtype='double', history=rstat )

    !雲密度
    call HistoryAddVariable(                         &
      & varname='DensCloud', dims=(/'x','z','t'/),   &
      & longname='DensCloud',                        &
      & units='kg/m^3',                              &
      & xtype='double', history=rstat )

    !飽和比
    call HistoryAddVariable(                         &
      & varname='SatRatio', dims=(/'x','z','t'/),    &
      & longname='SatRatio',                         &
      & units='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 ReStartFile_OutPut(        &
    &   Time, xz_PotTemp,               &
!   &   xz_PotTempSum,                  &
    &   xz_Exner, pz_VelX, xr_VelZ,     &
    &   xza_MixRt, xz_Km , xz_Kh,       &
    &   xz_DensCloud, xz_SatRatio       &
    & )
    !
    !リスタートファイルに予報変数を書き出す
    !
    
    !暗黙の型宣言禁止
    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_PotTempSum(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)  :: xz_DensCloud(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)  :: xz_SatRatio(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( 'PotTempDist', xz_PotTemp , rstat)
!    call HistoryPut( 'PotTemp', xz_PotTempSum , rstat)
    call HistoryPut( 'Km',      xz_Km      , rstat)
    call HistoryPut( 'Kh',      xz_Kh      , rstat)
    call HistoryPut( 'DensCloud',      xz_DensCloud      , rstat)
    call HistoryPut( 'SatRatio',      xz_SatRatio      , rstat)
    call HistoryPut( 'MixRt',   xza_MixRt  , rstat)    
    
  end subroutine ReStartFile_OutPut
  
  
  subroutine ReStartFile_Close
    !
    !リスタートファイルのクローズ
    !
    
    !モジュール読み込み
    use gt4_history
    
    !暗黙の型宣言禁止
    implicit none

    !ファイルを閉じる
    call HistoryClose(rstat)
    
  end subroutine ReStartFile_Close


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

    !暗黙の型宣言禁止
    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) :: xz_DensCloudN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_SatRatioN(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) :: xz_DensCloudB(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(out) :: xz_SatRatioB(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) = '^' // 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 = "DensCloud"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_DensCloudB = Var3D(:,:,1)
    xz_DensCloudN = Var3D(:,:,2)

    do t = 1, 2  
      name = "SatRatio"
      call HistoryGet( InitFile, name, Var3D(:,:,t), step(t) )    
    end do
    xz_SatRatioB = Var3D(:,:,1)
    xz_SatRatioN = 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
  
       
end module ReStartFileIO
