dcpam.f90

Path: main/dcpam.f90
Last Update: Mon Aug 21 17:01:09 JST 2006

Held and Suarez (1994) ベンチマークテストプログラム

Authors:Yasuhiro MORIKAWA
Version:$Id: dcpam.f90,v 1.3 2006/08/21 08:01:09 morikawa Exp $
Tag Name:$Name: dcpam3-20070608 $
Copyright:Copyright (C) GFD Dennou Club, 2005. All rights reserved.
License:See COPYRIGHT

Methods

dcpam  

Included Modules

nmlfile_mod constants_mod dycore_type_mod dycore_grid_mod dycore_in_mod dycore_out_mod dycore_dynamics_mod dycore_hs94forcing_mod dycore_time_mod dc_args dc_trace dc_message dc_string

Public Instance methods

Main Program :

dcpam のテストプログラムとしての Held and Suarez (1994) ベンチマークテスト を行なうプログラムです. 現在は main/init.f90 によって初期値を生成する 必要があります.

This procedure input/output NAMELIST#dcpam_nml .

[Source]

program dcpam
  !
  !dcpam のテストプログラムとしての Held and Suarez (1994) ベンチマークテスト
  !を行なうプログラムです. 現在は main/init.f90 によって初期値を生成する
  !必要があります.
  !

  !----- NAMELIST 名設定モジュール -----
  use nmlfile_mod, only: nmlfile_init, nmlfile_open, nmlfile_close

  !----- 物理定数設定モジュール -----
  use constants_mod, only: SecPerDay

  !----- 構造体参照モジュール -----
  use dycore_type_mod, only: DYCORE_VARS, DYCORE_DIMS, STRING, INTKIND, REKIND

  !----- 格子点取得モジュール -----
  use dycore_grid_mod, only: dycore_grid_init, dycore_grid_end, im, jm, km, nm

  !----- データI/Oモジュール -----
  use dycore_in_mod,   only: dycore_in_init, dycore_in_dims, dycore_in_vars, dycore_in_end
  use dycore_out_mod,  only: dycore_out_init, dycore_out_put, dycore_out_end

  !----- 力学過程モジュール (スペクトルモデル用) -----
  use dycore_dynamics_mod, only: dycore_dynamics_init, dycore_dynamics_prediction, dycore_dynamics_diagnostic, dycore_dynamics_diffusion, dycore_dynamics_end


  !----- Held and Suarez(1994) モジュール -----
  use dycore_hs94forcing_mod, only: dycore_hs94forcing_init, dycore_hs94forcing     , dycore_hs94forcing_end

  !----- 時間更新・タイムフィルターモジュール -----
  use dycore_time_mod, only: dycore_time_init, dycore_time_progress, dycore_time_end,  dycore_time_filter, StepInterval, OutputStep, CurrentLoop, CurrentTime, DelTime
!      ! &                     dycore_time_filter

  !----- コマンドライン引数解析ツール -----
  use dc_args, only: ARGS, Open, Debug, Help, Strict, Close

  !----- デバッグ・汎用ツール -----
  use dc_trace,  only: SetDebug, DbgMessage, BeginSub, EndSub
  use dc_message,only: MessageNotify
  use dc_string, only: toChar, Printf
  implicit none

  logical :: CalcDynamics    = .true. ! 力学コアの演算
  logical :: CalcDiffusion   = .true. ! 拡散項の演算
  logical :: CalcHS94forcing = .true. ! Held and Suarez(1994) による強制

  namelist /dcpam_nml/ CalcDynamics, CalcDiffusion, CalcHS94forcing
    ! 実験設定を行う。

  !-------------------------------------------------------------------
  !   変数定義
  !-------------------------------------------------------------------
  !----- 次元データ全てを包括 (時間次元のみ除く) -----
  type(DYCORE_DIMS)          :: Dims       ! 次元データ全種

  !----- 計算される格子点データ全てを包括 -----
  type(DYCORE_VARS)          :: Vars_b     ! 格子点データ全種(t-Δt)
  type(DYCORE_VARS)          :: Vars_n     ! 格子点データ全種(t)
  type(DYCORE_VARS)          :: Vars_a     ! 格子点データ全種(t+Δt)

  !----- 作業変数 -----
  !for time measure in temporarily
  real(REKIND) :: start_time_init, end_time_init, start_time_dynmcs, end_time_dynmcs, store_time_dynmcs = 0.0, start_time_diff, end_time_diff, store_time_diff = 0.0, start_time_diag, end_time_diag, store_time_diag = 0.0, start_time_hs94, end_time_hs94, store_time_hs94 = 0.0, total_time = 0.0 , remaining_time = 0.0
  integer(INTKIND) :: dateval(1:8)

  !for nmlfile_mod
  integer(INTKIND)   :: nmlstat, nmlunit
  logical            :: nmlreadable

  !for dc_args
  type(ARGS) :: arg

  character(*), parameter:: version = '$Name: dcpam3-20070608 $' // '$Id: dcpam.f90,v 1.3 2006/08/21 08:01:09 morikawa Exp $'
  character(STRING),  parameter:: subname = "dcpam"

continue
  !-------------------------------------------------------------------
  !   デバッグモード開始
  !-------------------------------------------------------------------
  call Open(arg)
  call Debug(arg)
   call Help(arg) 
   call Strict(arg)
  call Close(arg)

  call BeginSub(subname, version=version)

  !-------------------------------------------------------------------
  !   use CPU_TIME
  !-------------------------------------------------------------------
  call cpu_time(start_time_init)

  !-------------------------------------------------------------------
  !   NAMELIST file Setting
  !-------------------------------------------------------------------
  call nmlfile_init('dcpam.nml')

  !----------------------------------------------------------------
  !   Read init_nml
  !----------------------------------------------------------------
  call nmlfile_open(nmlunit, nmlreadable)
  if (nmlreadable) then
     read(nmlunit, nml=dcpam_nml, iostat=nmlstat)
     call DbgMessage('Stat of NAMELIST dcpam_nml Input is <%d>', i=(/nmlstat/))
     write(0, nml=dcpam_nml)
  else
     call DbgMessage('Not Read NAMELIST dcpam_nml')
     call MessageNotify('W', subname, 'Can not Read NAMELIST dcpam_nml. Force Use Default Value.')
  end if
  call nmlfile_close

  !-------------------------------------------------------------------
  !   格子点数の取得
  !-------------------------------------------------------------------
  call dycore_grid_init

  !-------------------------------------------------------------------
  !   初期値データ読み込み (gt4f90io の上位モジュール)
  !     (ファイル名や読み込む変数は NAMELIST で変更)
  !-------------------------------------------------------------------
  call dycore_in_init
  call dycore_in_dims( Dims          )  ! intent(out): 次元データ全種

  call dycore_in_vars( Vars_b      , Vars_n    )     ! intent(out): 格子点データ全種(t)

  !-------------------------------------------------------------------
  !   時間・ ステップに関する初期設定
  !-------------------------------------------------------------------
  call dycore_time_init

  !-------------------------------------------------------------------
  !   データ出力の初期設定
  !-------------------------------------------------------------------
  call dycore_out_init( Dims           ) ! intent(in): 次元データ全種

  !-------------------------------------------------------------------
  !   力学過程演算の初期設定
  !-------------------------------------------------------------------
  call dycore_dynamics_init( Dims              , Vars_a          )     ! intent(out): 格子点データ全種(t+Δt)

  !-------------------------------------------------------------------
  !   Held and Suarez(1994) モジュールの初期化
  !-------------------------------------------------------------------
  call dycore_hs94forcing_init( Dims            )     ! intent(in): 次元データ全種

  !-------------------------------------------------------------------
  !   use CPU_TIME
  !-------------------------------------------------------------------
  call cpu_time(end_time_init)
  call date_and_time(values=dateval)
  call Printf(6, 'Elapsed Time=<%r>, Start=<%r>, Stop=<%r>, ' // 'Now=<%d/%d/%d  %d:%d:%d>', r=(/end_time_init - start_time_init, start_time_init, end_time_init/), i=(/dateval(1), dateval(2), dateval(3), dateval(5), dateval(6), dateval(7)  /) )


  !-------------------------------------------------------------------
  !   (test) 初期値データがちゃんと読まれているかチェック
  !-------------------------------------------------------------------
!!$  call dycore_out_put    &
!!$       &  ( Vars_b   )  ! intent(in) : 格子点データ全種(t-Δt)
!!$  call dycore_out_put    &
!!$       &  ( Vars_n   )  ! intent(in) : 格子点データ全種(t)

  MainLoop : do
     call BeginSub( subname, 'MainLoop : CurrentLoop=<%d>, CurrentTime=<%f>', i=(/CurrentLoop/), d=(/CurrentTime/) )
     !----------------------------------------------------------------
     !   物理過程演算
     !----------------------------------------------------------------

     !----------------------------------------------------------------
     !   力学過程演算
     !----------------------------------------------------------------
     if (CalcDynamics) then
        call cpu_time(start_time_dynmcs)
        call dycore_dynamics_prediction( Dims          , Vars_b        , Vars_n        , Vars_a  )           ! intent(inout): 格子点データ全種(t+Δt)
        call cpu_time(end_time_dynmcs)
        store_time_dynmcs = store_time_dynmcs + (end_time_dynmcs - start_time_dynmcs)


        call cpu_time(start_time_diff)
        if (CalcDiffusion) then
           call dycore_dynamics_diffusion( Vars_b    , Vars_a  )      ! intent(inout): 格子点データ全種(t+Δt)
        else
           call DbgMessage('Not Consider Diffusion.')
        end if
        call cpu_time(end_time_diff)
        store_time_diff = store_time_diff + (end_time_diff - start_time_diff)


        call cpu_time(start_time_diag)
        !-----  診断出力  -----
        call dycore_dynamics_diagnostic( Dims            , Vars_a  )           ! intent(inout): 格子点データ全種(t+Δt)
        call cpu_time(end_time_diag)
        store_time_diag = store_time_diag + (end_time_diag - start_time_diag)

     else
        call cpu_time(start_time_dynmcs)
        call DbgMessage('Not Calculate Dynamics.')
        Vars_a%xyz_VelLon = Vars_b%xyz_VelLon  ! 速度経度成分
        Vars_a%xyz_VelLat = Vars_b%xyz_VelLat  ! 速度緯度成分
        Vars_a%xyz_Vor    = Vars_b%xyz_Vor   ! 渦度
        Vars_a%xyz_Div    = Vars_b%xyz_Div   ! 発散
        Vars_a%xyz_Temp   = Vars_b%xyz_Temp  ! 温度
        Vars_a%xyz_QVap   = Vars_b%xyz_QVap  ! 比湿
        Vars_a%xy_Ps      = Vars_b%xy_Ps     ! 地表面気圧
        call cpu_time(end_time_dynmcs)
        store_time_dynmcs = store_time_dynmcs + (end_time_dynmcs - start_time_dynmcs)


        call cpu_time(start_time_diff)
        if (CalcDiffusion) then
           call dycore_dynamics_diffusion( Vars_b    , Vars_a  )      ! intent(inout): 格子点データ全種(t+Δt)
        else
           call DbgMessage('Not Consider Diffusion.')
        end if
        call cpu_time(end_time_diff)
        store_time_diff = store_time_diff + (end_time_diff - start_time_diff)


        call cpu_time(start_time_diag)
        !-----  診断出力  -----
        call dycore_dynamics_diagnostic( Dims            , Vars_a  )           ! intent(inout): 格子点データ全種(t+Δt)
        call cpu_time(end_time_diag)
        store_time_diag = store_time_diag + (end_time_diag - start_time_diag)

     end if


     !----------------------------------------------------------------
     !   Held and Suarez(1994) の加熱、散逸の加算
     !----------------------------------------------------------------
     call cpu_time(start_time_hs94)
     if (CalcHS94forcing) then
        call dycore_hs94forcing( Vars_b     , Vars_n     , Vars_a   )     ! intent(inout): 格子点データ全種(t+Δt)
     else
        call DbgMessage('Not Calculate Held-Suarez Forcing.')
     end if
     call cpu_time(end_time_hs94)
     store_time_hs94 = store_time_hs94 + (end_time_hs94 - start_time_hs94)

     !----------------------------------------------------------------
     !   タイムフィルター
     !----------------------------------------------------------------
     call dycore_time_filter( Vars_b           , Vars_n           , Vars_a  )            ! intent(in)   : 格子点データ全種(t+Δt)

     !----------------------------------------------------------------
     !   データの出力
     !----------------------------------------------------------------
     call dycore_out_put( Vars_n  )            ! intent(in): 格子点データ全種(t)

!!$     !----------------------------------------------------------------
!!$     !   リスタートファイル出力
!!$     !----------------------------------------------------------------
!!$     call dycore_restart_put( &
!!$       & Vars_n           , & ! intent(in): 格子点データ全種(t)
!!$       & Vars_a  )            ! intent(in): 格子点データ全種(t+Δt)

     !----------------------------------------------------------------
     !   時刻更新
     !----------------------------------------------------------------
     call dycore_time_progress( Vars_b           , Vars_n           , Vars_a  )            ! intent(inout): 格子点データ全種(t+Δt)

     !-------------------------------------------------------------------
     !   use CPU_TIME
     !-------------------------------------------------------------------
     if ( mod(CurrentLoop, StepInterval) == 0 ) then
        call Printf(6, '')
        call Printf(6, '-- ***              ***  ')
        call Printf(6, '-- *** TIME         ***  ')
        call Printf(6, '-- ***      KEEPING ***  ')
        call Printf(6, '-- ***              ***  ')
        call Printf(6, '--CurrentLoop=<%d>, CurrentTime=<%f>, Day=<%f>', i=(/CurrentLoop/), d=(/CurrentTime, CurrentTime/SecPerDay /) )

        call Printf(6, '--Remaining Loop  ' // ': RemainingLoop=<%d>, RemainingTime=<%f>, Day=<%f>', i=(/StepInterval * OutputStep - CurrentLoop/), d=(/StepInterval * OutputStep * DelTime - CurrentTime , ( StepInterval * OutputStep * DelTime - CurrentTime ) / SecPerDay /)    )

        call cpu_time(total_time)
        remaining_time = ( total_time / real(CurrentLoop, REKIND) ) * real( StepInterval * OutputStep - CurrentLoop, REKIND )

        call Printf(6, '--Remaining Real Time Forecast  ' // ': RemainingRealTime=<%r>, Min=<%r>, Hour=<%r>, Day=<%r>', r=(/remaining_time, remaining_time / 6.0e1, remaining_time/ 3.6e3, remaining_time/ 8.64e4 /) )

        call date_and_time(values=dateval)
        call Printf(6, '--Current Time <%d/%d/%d  %d:%d:%d>', i=(/dateval(1), dateval(2), dateval(3), dateval(5), dateval(6), dateval(7)  /) )

        call Printf(6, '--[Categoly]  : [Stored]       ' // '[Current]      [Start]        [Stop] '  )

        call Printf(6, '--Dynamics    : <%r> <%r> <%r> <%r> ', r=(/store_time_dynmcs, end_time_dynmcs - start_time_dynmcs, start_time_dynmcs, end_time_dynmcs/)   )

        call Printf(6, '--Diffusion   : <%r> <%r> <%r> <%r> ', r=(/store_time_diff, end_time_diff - start_time_diff, start_time_diff, end_time_diff/)   )

        call Printf(6, '--Diagnostic  : <%r> <%r> <%r> <%r> ', r=(/store_time_diag, end_time_diag - start_time_diag, start_time_diag, end_time_diag/)   )

        call Printf(6, '--hs94forcing : <%r> <%r> <%r> <%r> ', r=(/store_time_hs94, end_time_hs94 - start_time_hs94, start_time_hs94, end_time_hs94/)   )
     end if

     !----------------------------------------------------------------
     !   ループから抜け出すかどうかの判定
     !----------------------------------------------------------------
     if ( CurrentLoop > StepInterval * OutputStep ) then
        exit MainLoop
     endif

     call EndSub(subname)
  enddo MainLoop

  call dycore_out_end

  call EndSub(subname)
end program dcpam

[Validate]