!= 力学過程 (スペクトル法, Arakawa and Suarez (1983))
!
!= Dynamical Core (Spectral method, Arakawa and Suarez (1983))
!
! Authors:: Yasuhiro MORIKAWA
! Version:: $Id: dyn_spectral_as83.f90,v 1.7 2007/07/31 06:53:16 morikawa Exp $
! Tag Name:: $Name: dcpam4-20070731-1 $
! Copyright:: Copyright (C) GFD Dennou Club, 2007. All rights reserved.
! License:: See COPYRIGHT[link:../../COPYRIGHT]
!
module dyn_spectral_as83
!
!= 力学過程 (スペクトル法, Arakawa and Suarez (1983))
!
!= Dynamical Core (Spectral method, Arakawa and Suarez (1983))
!
! Note that Japanese and English are described in parallel.
!
! 力学過程を演算するモジュールです.
! 水平離散化にスペクトル法を, 鉛直離散化には Arakawa and Suarez (1983)
! を用いています.
!
! This is a dynamical core module. Spectral method (for horizontal) and
! Arakawa and Suarez (1983) method (for vertical) are used.
!
!== Procedures List
!
! Create :: DYNSPAS83 型変数の初期化
! GetAxes :: 座標値の取得
! Dynamics :: 力学過程の演算
! Close :: DYNSPAS83 型変数の終了処理
! PutLine :: DYNSPAS83 型変数に格納されている情報の印字
! initialized :: DYNSPAS83 型変数が初期設定されているか否か
! ------------ :: ------------
! Create :: Constructor of "DYNSPAS83"
! GetAxes :: Get data of axes
! Dynamics :: Calculation of dynamical core
! Close :: Deconstructor of "DYNSPAS83"
! PutLine :: Print information of "DYNSPAS83"
! initialized :: Check initialization of "DYNSPAS83"
!
!== Usage
!
! 主プログラムの始めに, DYNSPAS83 型の変数を Create で初期化します.
! 次に, 時間積分ループ内で渦度, 発散, 温度, 比湿, 地表面気圧
! ( $ t $ および $ t-\Delta t $ ) を Dynamics に与えてください.
! それらの値を用いて力学過程の演算が行われ,
! $ t+\Delta t $ の東西風速, 南北風速, 渦度, 発散, 温度, 比湿,
! 地表面気圧が得られます.
! 最後に, DYNSPAS83 型の変数の終了処理を Close にて行います.
! 主プログラムで座標値が必要となる場合には, GetAxes を用いて
! 座標値を取得してください.
!
! In the begging of program, initialize "DYNSPAS83" by Create.
! Next, in the time integration loop,
! provide vorticity, divergence, temperature, specific humidity,
! surface pressure on $ t $ and $ t-\Delta t $ to Dynamics.
! Then the physical values at $ t+\Delta t $ is returned.
! Finally, terminate "DYNSPAS83" by Close.
! If you need data of axes, get that by GetAxes.
!
!== References
!
! * Arakawa, A., Suarez, M. J., 1983:
! Vertical differencing of the primitive equations
! in sigma coordinates.
! Mon. Wea. Rev., 111, 34--35.
use dyn_spectral, only: DYNSP
use dyn_as83, only: DYNAS83
implicit none
private
public:: DYNSPAS83, Create, Close, PutLine, initialized
public:: GetAxes, Dynamics, VorDiv2UV, UV2VorDiv
type DYNSPAS83
!
! まず, Create で "DYNSPAS83" 型の変数を初期設定して下さい.
! 初期化された "DYNSPAS83" 型の変数を再度利用する際には,
! Close によって終了処理を行ってください.
!
! Initialize "DYNSPAS83" variable by "Create" before usage.
! If you reuse "DYNSPAS83" variable again for another application,
! terminate by "Close".
!
logical:: initialized = .false. ! 初期設定フラグ.
! Initialization flag
integer:: nmax ! 最大全波数.
! Maximum truncated wavenumber
integer:: imax ! 経度格子点数.
! Number of grid points in longitude
integer:: jmax ! 緯度格子点数.
! Number of grid points in latitude
integer:: kmax ! 鉛直層数.
! Number of vertical level
type(DYNSP):: dyn_sp
type(DYNAS83):: dyn_as
end type DYNSPAS83
character(*), parameter:: version = &
& '$Name: dcpam4-20070731-1 $' // &
& '$Id: dyn_spectral_as83.f90,v 1.7 2007/07/31 06:53:16 morikawa Exp $'
interface Create
module procedure DynSpectralAS83Create
end interface
interface Close
module procedure DynSpectralAS83Close
end interface
interface PutLine
module procedure DynSpectralAS83PutLine
end interface
interface initialized
module procedure DynSpectralAS83Initialized
end interface
interface NmlRead
module procedure DynSpectralAS83NmlRead
end interface
interface GetAxes
module procedure DynSpectralAS83GetAxes
end interface
interface Dynamics
module procedure DynSpectralAS83Dynamics
end interface
interface VorDiv2UV
module procedure DynSpectralAS83VorDiv2UV
end interface
interface UV2VorDiv
module procedure DynSpectralAS83UV2VorDiv
end interface
contains
subroutine DynSpectralAS83Create( dyn_sp_as, &
& nmax, imax, jmax, kmax, &
& PI, RPlanet, Grav, Omega, Cp, RAir, EpsVT, VisOrder, EFoldTime, &
& DelTime, &
& xy_Phis, &
& openmp_threads, &
& r_SigmaSet, &
& nmlfile, err)
!
! 引数 *dyn_sp_as* に力学過程の設定を行う.
! 他のサブルーチンを使用する前に必ずこのサブルーチンによって
! DYNSPAS83 型の変数を初期化してください.
!
! Configure the settings for dynamical core to *dyn_sp_as*.
! Before other subroutines are used, initialize "DYNSPAS83"
! variable.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, DP, STDOUT
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_EALREADYINIT, &
& DCPAM_ENEGATIVE
use dyn_spectral, only: Create, &
& GetCoriolis, GetSpectralCoeff, GetDiffCoeff, GetPhis
use dyn_as83, only: Create
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
integer, intent(in):: nmax ! 最大全波数.
! Maximum truncated wavenumber
integer, intent(in):: imax ! 経度格子点数.
! Number of grid points in longitude
integer, intent(in):: jmax ! 緯度格子点数.
! Number of grid points in latitude
integer, intent(in):: kmax ! 鉛直層数.
! Number of vertical level
real(DP), intent(in):: PI ! $ \pi $ . 円周率. Circular constant
real(DP), intent(in):: RPlanet ! $ a $ . 惑星半径. Radius of planet
real(DP), intent(in):: Omega ! $ \Omega $ . 回転角速度. Angular velocity
real(DP), intent(in):: Grav ! $ g $ . 重力加速度. Gravitational acceleration
real(DP), intent(in):: Cp ! $ C_p $ . 大気定圧比熱. Specific heat of air at constant pressure
real(DP), intent(in):: RAir ! $ R $ . 大気気体定数. Gas constant of air
real(DP), intent(in):: EpsVT ! $ 1/\epsilon_v - 1 $ .
integer, intent(in):: VisOrder ! 超粘性の次数. Order of hyper-viscosity
real(DP), intent(in):: EFoldTime ! 最大波数に対する e-folding time. E-folding time for maximum wavenumber
real(DP), intent(in):: DelTime ! $ \Delta t $ . タイムステップ. Time step
real(DP), intent(in), optional:: xy_Phis (0:imax-1, 0:jmax-1)
! $ \Phi_s $ . 地表ジオポテンシャル.
! Surface geo-potential
integer, intent(in), optional:: openmp_threads
! OPENMP での最大スレッド数.
! openmp_threads に 1 より大きな値を指定すれば
! ISPACK[http://www.gfd-dennou.org/library/ispack/]
! の球面調和函数変換 OPENMP 並列計算
! ルーチンが用いられる. 並列計算を実行するには,
! 実行時に環境変数 OMP_NUM_THREADS
! を openmp_threads 以下の数字に設定する
! 等のシステムに応じた準備が必要となる.
!
! openmp_threads に 1 より大きな値を
! 指定しなければ並列計算ルーチンは呼ばれない.
!
! Maximum number of threads in OPENMP.
! If integer more than 1 is given,
! spherical harmonics transration subroutines
! with OPENMP parallel computation in
! ISPACK[http://www.gfd-dennou.org/library/ispack/]
! is used. In practice, specify environment
! variable (OMP_NUM_THREADS etc.) for parallel
! computation.
!
! If number less than 2 is given,
! subroutines of parallel computation
! is not called.
real(DP), intent(in), optional:: r_SigmaSet (:)
! デフォルトでは, kmax の値に応じ,
! 自動的に $ \sigma $ レベル (半整数)
! は設定される (kmax がある特定の値のみ).
! $ \sigma $ レベル (半整数).
! を明示的に設定する必要がある場合,
! この引数を与える.
!
! By default, half $ \sigma $ level is
! specified automatically according to
! the value of kmax (only the value is
! same as particular values).
! If the half $ \sigma $ level is specified
! manually, give this argument.
character(*), intent(in), optional :: nmlfile
! NAMELIST ファイルの名称.
! この引数に空文字以外を与えた場合,
! 指定されたファイルから
! NAMELIST 変数群を読み込みます.
! ファイルを読み込めない場合にはエラーを
! 生じます.
!
! NAMELIST file name.
! If nonnull character is specified to
! this argument,
! NAMELIST group name is loaded from the
! file.
! If the file can not be read,
! an error occurs.
!
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
!-----------------------------------
! 作業変数
! Work variables
real(DP):: xy_Cori (0:imax-1, 0:jmax-1)
! $ f\equiv 2\Omega\sin\varphi $ .
! コリオリパラメータ. Coriolis parameter
integer:: nmo (1:2, 0:nmax, 0:nmax)
! スペクトルの添字順番.
! Spectral subscript expression
real(DP):: wz_rn ((nmax+1)**2, 0:kmax-1)
! $ -n \times (n+1) $ . ラプラシアンの係数.
! Laplacian coefficient
real(DP):: wz_DiffVorDiv ((nmax+1)**2, 0:kmax-1)
! $ -K_{HD} [(-1)^{N_D/2}\nabla^{N_D}- (2/a^2)^{N_D/2}] $ .
! 運動量水平拡散係数.
! Coefficient of horizontal momentum diffusion
real(DP):: wz_DiffTherm ((nmax+1)**2, 0:kmax-1)
! $ -(-1)^{N_D/2}K_{HD}\nabla^{N_D} $ .
! 熱, 水水平拡散係数.
! Coefficient of horizontal thermal and water diffusion
real(DP):: w_Phis ((nmax+1)**2)
! $ \Phi_s $ . 地表ジオポテンシャル.
! Surface geo-potential
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83Create'
continue
call BeginSub(subname, version)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (dyn_sp_as % initialized) then
stat = DCPAM_EALREADYINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! 引数の正当性のチェック
! Validation of arguments
!-----------------------------------------------------------------
if (nmax < 1) then
stat = DCPAM_ENEGATIVE
cause_c = 'nmax'
goto 999
end if
if (imax < 1) then
stat = DCPAM_ENEGATIVE
cause_c = 'imax'
goto 999
end if
if (jmax < 1) then
stat = DCPAM_ENEGATIVE
cause_c = 'jmax'
goto 999
end if
if (kmax < 1) then
stat = DCPAM_ENEGATIVE
cause_c = 'kmax'
goto 999
end if
!-------------------------
! NAMELIST からの値
! Values from NAMELIST
!!$ if ( present_and_not_empty(nmlfile) ) then
!!$ call MessageNotify( 'M', subname, &
!!$ & 'Loading NAMELIST file "%c" ...', &
!!$ & c1=trim(nmlfile) )
!!$ call NmlRead ( nmlfile = nmlfile, & ! (in)
!!$ & param_i = dyn_sp_as % param_i, & ! (inout)
!!$ & param_r = dyn_sp_as % param_r, & ! (inout)
!!$ & param_c_ = dyn_sp_as % param_c, & ! (inout)
!!$ & err = err ) ! (out)
!!$ if ( present_and_true(err) ) then
!!$ call MessageNotify( 'W', subname, &
!!$ & '"%c" can not be read.', &
!!$ & c1=trim(nmlfile) )
!!$ stat = DC_ENOFILEREAD
!!$ cause_c = nmlfile
!!$ goto 999
!!$ end if
!!$ end if
!-----------------------------------------------------------------
! 波数・格子点の設定
! Configure wave number and grid point
!-----------------------------------------------------------------
dyn_sp_as % nmax = nmax
dyn_sp_as % imax = imax
dyn_sp_as % jmax = jmax
dyn_sp_as % kmax = kmax
!-----------------------------------------------------------------
! スペクトル法の設定
! Configure the settings for spectral method
!-----------------------------------------------------------------
call Create( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& nmax = nmax, imax = imax, jmax = jmax, kmax = kmax, & ! (in)
& PI = PI, RPlanet = RPlanet, Grav = Grav, Omega = Omega, & ! (in)
& Cp = Cp, RAir = RAir, VisOrder = VisOrder, & ! (in)
& EFoldTime = EFoldTime, & ! (in)
& DelTime = DelTime, & ! (in)
& xy_Phis = xy_Phis, & ! (in)
& openmp_threads = openmp_threads, & ! (in)
& nmlfile = nmlfile ) ! (in)
!-----------------------------------------------------------------
! Arakawa and Suarez (1983) の設定
! Configure the settings for Arakawa and Suarez (1983)
!-----------------------------------------------------------------
!---------------------------------------------
! コリオリパラメータの取得
! Get Coriolis parameter
call GetCoriolis( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xy_Cori = xy_Cori ) ! (out)
!---------------------------------------------
! スペクトル係数の取得
! Get spectral coefficient
call GetSpectralCoeff( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& nmo = nmo, wz_rn = wz_rn ) ! (out)
!---------------------------------------------
! 拡散係数の取得
! Get diffusion coefficient
call GetDiffCoeff( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& wz_DiffVorDiv = wz_DiffVorDiv, wz_DiffTherm = wz_DiffTherm ) ! (out)
!---------------------------------------------
! 地形データ (地表 $ \Phi $ )の取得
! Get geography data (surface $ \Phi $ )
call GetPhis( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& w_Phis = w_Phis ) ! (out)
call Create( dyn_as = dyn_sp_as % dyn_as, & ! (inout)
& nmax = nmax, imax = imax, jmax = jmax, kmax = kmax, & ! (in)
& PI = PI, RPlanet = RPlanet, Grav = Grav, Omega = Omega, & ! (in)
& Cp = Cp, RAir = RAir, EpsVT = EpsVT, & ! (in)
& EFoldTime = EFoldTime, DelTime = DelTime, & ! (in)
& xy_Cori = xy_Cori, & ! (in)
& nmo = nmo, wz_rn = wz_rn, & ! (in)
& wz_DiffVorDiv = wz_DiffVorDiv, & ! (in)
& wz_DiffTherm = wz_DiffTherm, & ! (in)
& r_SigmaSet = r_SigmaSet, w_Phis = w_Phis, & ! (in)
& nmlfile = nmlfile ) ! (in)
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
dyn_sp_as % initialized = .true.
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83Create
subroutine DynSpectralAS83Close( dyn_sp_as, err )
!
! DYNSPAS83 型の変数の終了処理を行います.
! なお, 与えられた *dyn_sp_as* が Create によって初期設定
! されていない場合, プログラムはエラーを発生させます.
!
! Deconstructor of "DYNSPAS83".
! Note that if *dyn_sp_as* is not initialized by Create yet,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dyn_spectral, only: Close
use dyn_as83, only: Close
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83Close'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (.not. dyn_sp_as % initialized) then
stat = DCPAM_ENOTINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! スペクトル法と Arakawa and Suarez (1983) の設定の消去
! Clear the settings for spectral method and Arakawa and Suarez (1983)
!-----------------------------------------------------------------
call Close(dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& err = err) ! (out)
call Close(dyn_as = dyn_sp_as % dyn_as, & ! (inout)
& err = err) ! (out)
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
dyn_sp_as % initialized = .false.
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83Close
subroutine DynSpectralAS83PutLine( dyn_sp_as, unit, indent, err )
!
! 引数 *dyn_sp_as* に設定されている情報を印字します.
! デフォルトではメッセージは標準出力に出力されます.
! *unit* に装置番号を指定することで, 出力先を変更することが可能です.
!
! Print information of *dyn_sp_as*.
! By default messages are output to standard output.
! Unit number for output can be changed by *unit* argument.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dc_string, only: Printf
use dyn_spectral, only: PutLine
use dyn_as83, only: PutLine
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
integer, intent(in), optional:: unit
! 出力先の装置番号.
! デフォルトの出力先は標準出力.
!
! Unit number for output.
! Default value is standard output.
character(*), intent(in), optional:: indent
! 表示されるメッセージの字下げ.
!
! Indent of displayed messages.
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
integer:: stat
character(STRING):: cause_c
integer:: out_unit
integer:: indent_len
character(STRING):: indent_str
character(*), parameter:: subname = 'DynSpectralAS83PutLine'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (present(unit)) then
out_unit = unit
else
out_unit = STDOUT
end if
indent_len = 0
indent_str = ''
if (present(indent)) then
if (len(indent) /= 0) then
indent_len = len(indent)
indent_str(1:indent_len) = indent
end if
end if
!-----------------------------------------------------------------
! "DYNSPAS83" の設定の印字
! Print the settings for "DYNSPAS83"
!-----------------------------------------------------------------
if (dyn_sp_as % initialized) then
call Printf(out_unit, &
& indent_str(1:indent_len) // &
& '#')
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83PutLine
logical function DynSpectralAS83Initialized( dyn_sp_as ) result(result)
!
! *dyn_sp_as* が初期設定されている場合には .true. が,
! 初期設定されていない場合には .false. が返ります.
!
! If *dyn_sp_as* is initialized, .true. is returned.
! If *dyn_sp_as* is not initialized, .false. is returned.
!
implicit none
type(DYNSPAS83), intent(in):: dyn_sp_as
continue
result = dyn_sp_as % initialized
end function DynSpectralAS83Initialized
subroutine DynSpectralAS83NmlRead( nmlfile, &
!!$ & param_i, param_r, param_c_, &
& err )
!
! NAMELIST ファイル *nmlfile* から値を入力するための
! 内部サブルーチンです. Create 内で呼び出されることを
! 想定しています.
!
! 値が NAMELIST ファイル内で指定されていない場合には,
! 入力された値がそのまま返ります.
!
! なお, *nmlfile* に空文字が与えられた場合, または
! 与えられた *nmlfile* を読み込むことができない場合,
! プログラムはエラーを発生させます.
!
! This is an internal subroutine to input values from
! NAMELIST file *nmlfile*. This subroutine is expected to be
! called by "Create".
!
! A value not specified in NAMELIST file is returned
! without change.
!
! If *nmlfile* is empty, or *nmlfile* can not be read,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: DP, STRING, TOKEN, STDOUT
use dc_iounit, only: FileOpen
use dc_message, only: MessageNotify
use dc_present, only: present_and_true
use dc_error, only: StoreError, DC_NOERR, DC_ENOFILEREAD
implicit none
character(*), intent(in):: nmlfile
! NAMELIST ファイルの名称.
! NAMELIST file name
!!$ integer, intent(inout):: param_i
!!$ real(DP), intent(inout):: param_r
!!$ character(*), intent(inout):: param_c_
!!$ character(TOKEN):: param_c
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
!!$ namelist /dyn_spectral_as83_nml/ &
!!$ & param_i, param_r, param_c
! dyn_spectral_as83 モジュール用
! NAMELIST 変数群名.
!
! NAMELIST group name for
! 'dyn_spectral_as83' module.
!-----------------------------------
! 作業変数
! Work variables
integer:: stat
character(STRING):: cause_c
integer:: unit_nml ! NAMELIST ファイルオープン用装置番号.
! Unit number for NAMELIST file open
integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT.
! IOSTAT of NAMELIST read
character(*), parameter:: subname = 'DynSpectralAS83NmlRead'
continue
call BeginSub( subname )
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 文字型引数を NAMELIST 変数群へ代入
! Substitute character arguments to NAMELIST group
!-----------------------------------------------------------------
!!$ param_c = param_c_
!----------------------------------------------------------------
! NAMELIST ファイルのオープン
! Open NAMELIST file
!----------------------------------------------------------------
call FileOpen( unit = unit_nml, & ! (out)
& file = nmlfile, mode = 'r', & ! (in)
& err = err ) ! (out)
if ( present_and_true(err) ) then
stat = DC_ENOFILEREAD
cause_c = nmlfile
goto 999
end if
!-----------------------------------------------------------------
! NAMELIST 変数群の取得
! Get NAMELIST group
!-----------------------------------------------------------------
!!$ read( unit = unit_nml, & ! (in)
!!$ & nml = dyn_spectral_as83_nml, iostat = iostat_nml ) ! (out)
!!$ if ( iostat_nml == 0 ) then
!!$ call MessageNotify( 'M', subname, &
!!$ & 'NAMELIST group "%c" is loaded from "%c".', &
!!$ & c1='dyn_spectral_as83_nml', c2=trim(nmlfile) )
!!$ write(STDOUT, nml = dyn_spectral_as83_nml)
!!$ else
!!$ call MessageNotify( 'W', subname, &
!!$ & 'NAMELIST group "%c" is not found in "%c" (iostat=%d).', &
!!$ & c1='dyn_spectral_as83_nml', c2=trim(nmlfile), &
!!$ & i=(/iostat_nml/) )
!!$ end if
close( unit_nml )
!-----------------------------------------------------------------
! NAMELIST 変数群を文字型引数へ代入
! Substitute NAMELIST group to character arguments
!-----------------------------------------------------------------
!!$ param_c_ = param_c
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError( stat, subname, err, cause_c )
call EndSub( subname )
end subroutine DynSpectralAS83NmlRead
subroutine DynSpectralAS83GetAxes( dyn_sp_as, &
& x_Lon, y_Lat, z_Sigma, r_Sigma, z_DelSigma, &
& err )
!
! *dyn_sp_as* に格納されている座標値を返します.
!
! なお, 与えられた *dyn_sp_as* が Create によって初期設定
! されていない場合, プログラムはエラーを発生させます.
!
! Return data of axes stored in *dyn_sp_as*.
!
! If *dyn_sp_as* is not initialized by Create yet,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT, DP
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dyn_spectral, only: GetAxes
use dyn_as83, only: GetAxes
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
real(DP), intent(out):: x_Lon (0:dyn_sp_as%imax-1)
! 経度. Longitude
real(DP), intent(out):: y_Lat (0:dyn_sp_as%jmax-1)
! 緯度. Latitude
real(DP), intent(out):: z_Sigma (0:dyn_sp_as%kmax-1)
! $ \sigma $ レベル (整数).
! Full $ \sigma $ level
real(DP), intent(out):: r_Sigma (0:dyn_sp_as%kmax)
! $ \sigma $ レベル (半整数).
! Half $ \sigma $ level
real(DP), intent(out), optional:: z_DelSigma (0:dyn_sp_as%kmax-1)
! $ \Delta \sigma $ (整数).
! $ \Delta \sigma $ (Full)
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83Sample'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (.not. dyn_sp_as % initialized) then
stat = DCPAM_ENOTINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! 緯度経度座標値の取得
! Get data of latitude and longitude
!-----------------------------------------------------------------
call GetAxes( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& x_Lon = x_Lon, y_Lat = y_Lat ) ! (out)
!-----------------------------------------------------------------
! 鉛直レベル座標値の取得
! Get data of vertical level
!-----------------------------------------------------------------
call GetAxes( dyn_as = dyn_sp_as % dyn_as, & ! (inout)
& z_Sigma = z_Sigma, r_Sigma = r_Sigma, & ! (out)
& z_DelSigma = z_DelSigma ) ! (out)
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83GetAxes
subroutine DynSpectralAS83Dynamics( dyn_sp_as, &
& xyz_VorB, xyz_DivB, xyz_TempB, xyz_QVapB, xy_PsB, &
& xyz_VorN, xyz_DivN, xyz_TempN, xyz_QVapN, xy_PsN, &
& xyz_VorA, xyz_DivA, xyz_TempA, xyz_QVapA, xy_PsA, &
& err )
!
! 力学過程の演算を行い, 与えられた $ t-\Delta t $ および $ t $ の
! 渦度, 発散, 温度, 比湿, 地表面気圧から,
! $ t+\Delta t $ の
! 渦度, 発散, 温度, 比湿, 地表面気圧 を返します.
!
! 時間積分法にはリープフロッグスキームを用いています.
! さらに $ \Delta t $ を大きくとるために, 重力波項に
! セミインプリシット法を適用しています.
!
! なお, 与えられた *dyn_sp_as* が Create によって初期設定
! されていない場合, プログラムはエラーを発生させます.
!
! Calculating dynamical core,
! from vorticity, divergence, temperature, specific humidity,
! surface pressure at $ t-\Delta t $ and $ t $, the
! physical values at $ t+\Delta t $ are returned.
!
! Leap-frog scheme is used as time integration.
! In addition, semi-implicit scheme is applied
! for extension of $ \Delta t $ .
!
! If *dyn_sp_as* is not initialized by Create yet,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT, DP
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dyn_spectral, only: GradPi, VorDiv2UV, Tendency, &
& Spectral2Grid, Grid2Spectral, &
& DiffusionCorrectTemp, DiffusionVorDiv
use dyn_as83, only: NonLinearOnGrid, TimeIntegration
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
real(DP), intent(in):: xyz_VorB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ \zeta (t-\Delta t) $ . 渦度. Vorticity
real(DP), intent(in):: xyz_DivB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ D (t-\Delta t) $ . 発散. Divergence
real(DP), intent(in):: xyz_TempB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ T (t-\Delta t) $ . 温度. Temperature
real(DP), intent(in):: xy_PsB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ P_s (t-\Delta t) $ . 地表面気圧. Surface pressure
real(DP), intent(in):: xyz_QVapB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ q (t-\Delta t) $ . 比湿. Specific humidity
real(DP), intent(in):: xyz_VorN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ \zeta (t) $ . 渦度. Vorticity
real(DP), intent(in):: xyz_DivN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ D (t) $ . 発散. Divergence
real(DP), intent(in):: xyz_TempN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ T (t) $ . 温度. Temperature
real(DP), intent(in):: xy_PsN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ P_s (t) $ . 地表面気圧. Surface pressure
real(DP), intent(in):: xyz_QVapN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ q (t) $ . 比湿. Specific humidity
real(DP), intent(out):: xyz_VorA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ \zeta (t+\Delta t) $ . 渦度. Vorticity
real(DP), intent(out):: xyz_DivA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ D (t+\Delta t) $ . 発散. Divergence
real(DP), intent(out):: xyz_TempA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ T (t+\Delta t) $ . 温度. Temperature
real(DP), intent(out):: xy_PsA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ P_s (t+\Delta t) $ . 地表面気圧. Surface pressure
real(DP), intent(out):: xyz_QVapA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ q (t+\Delta t) $ . 比湿. Specific humidity
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
!-----------------------------------
! 作業変数
! Work variables
real(DP):: xy_GradLonPiN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ \DD{\pi}{x} (t) $
real(DP):: xy_GradLatPiN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ \DD{\pi}{y} (t) $
real(DP):: xyz_UN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ U (t) $ . 東西風速. Zonal wind.
real(DP):: xyz_VN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ V (t) $ . 南北風速. Meridional wind.
real(DP):: xyz_UA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ U (t+\Delta t) $ . 東西風速. Zonal wind
real(DP):: xyz_VA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ V (t+\Delta t) $ . 南北風速. Meridional wind
real(DP):: xyz_UAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ UA (t) $ . 東西運動量移流項.
! Zonal advection of momentum
real(DP):: xyz_VAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ VA (t) $ . 南北運動量移流項.
! Meridional advection of momentum
real(DP):: xyz_DTempDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ H (t) $ . 温度時間変化項.
! Temperature tendency
real(DP):: xyz_DQVapDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ R (t) $ . 比湿時間変化項.
! Specific humidity tendency
real(DP):: xyz_KEN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ E (t) $ . 運動エネルギー項.
! Kinematic energy
real(DP):: xyz_TempUAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ UT (t) $ . 温度東西移流項.
! Zonal advection of temperature
real(DP):: xyz_TempVAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ VT (t) $ . 温度南北移流項.
! Meridional advection of temperature
real(DP):: xy_DPiDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1)
! $ Z $ . 地表面気圧時間変化項.
! Surface pressure tendency
real(DP):: xyz_QVapUAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ Uq (t) $ . 比湿東西移流項.
! Zonal advection of specific humidity
real(DP):: xyz_QVapVAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! $ Vq (t) $ . 比湿南北移流項.
! Meridional advection of specific humidity
real(DP):: wz_DVorDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \DD{\zeta}{t} (t) $ . 渦度変化 (スペクトル).
! Vorticity tendency (spectral)
real(DP):: wz_DDivDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \DD{D}{t} (t) $ . 発散変化 (スペクトル).
! Divergence tendency (spectral)
real(DP):: wz_DTempDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \DD{T}{t} (t) $ . 温度変化 (スペクトル).
! Temperature tendency (spectral)
real(DP):: wz_DQVapDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \DD{q}{t} (t) $ . 比湿変化 (スペクトル).
! Specific humidity tendency (spectral)
real(DP):: w_DPiDtN ((dyn_sp_as%nmax+1)**2)
! $ \DD{Ps}{t} (t) $ . 地表面気圧変化 (スペクトル).
! Surface pressure tendency (spectral)
real(DP):: wz_VorB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \zeta (t-\Delta t) $ . 渦度 (スペクトル).
! Vorticity (spectral)
real(DP):: wz_DivB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ D (t-\Delta t) $ . 発散 (スペクトル).
! Divergence (spectral)
real(DP):: wz_TempB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ T (t-\Delta t) $ . 温度 (スペクトル).
! Temperature (spectral)
real(DP):: w_PiB ((dyn_sp_as%nmax+1)**2)
! $ \pi = \ln P_s (t-\Delta t) $ . 地表面気圧 (スペクトル).
! Surface pressure (spectral)
real(DP):: wz_QVapB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ q (t-\Delta t) $ . 比湿 (スペクトル).
! Specific humidity (spectral)
real(DP):: wz_VorA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \zeta (t+\Delta t) $ . 渦度 (スペクトル).
! Vorticity (spectral)
real(DP):: wz_DivA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ D (t+\Delta t) $ . 発散 (スペクトル).
! Divergence (spectral)
real(DP):: wz_TempA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ T (t+\Delta t) $ . 温度 (スペクトル).
! Temperature (spectral)
real(DP):: w_PiA ((dyn_sp_as%nmax+1)**2)
! $ \pi = \ln P_s (t+\Delta t) $ . 地表面気圧 (スペクトル).
! Surface pressure (spectral)
real(DP):: wz_QVapA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ q (t+\Delta t) $ . 比湿 (スペクトル).
! Specific humidity (spectral)
real(DP):: wz_VorDiffA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \mathscr{D}(\zeta) $ .
! 運動量水平拡散による渦度変化 (スペクトル).
! Vorticity tendency by
! horizontal momentum diffusion (spectral)
real(DP):: wz_DivDiffA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1)
! $ \mathscr{D}(D) $ .
! 運動量水平拡散による発散変化 (スペクトル).
! Divergence tendency by
! horizontal momentum diffusion (spectral)
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83Dynamics'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (.not. dyn_sp_as % initialized) then
stat = DCPAM_ENOTINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! 時間変化項の計算
! Calculate tendency terms
!-----------------------------------------------------------------
!---------------------------------------------
! 地表面気圧の空間変化項の計算
! Calculate spatial surface pressure tendency
!
call GradPi( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xy_Ps = xy_PsN, & ! (in)
& xy_GradLonPi = xy_GradLonPiN, &
& xy_GradLatPi = xy_GradLatPiN ) ! (out)
!---------------------------------------------
! 渦度発散から風速の計算
! Calculate wind velocity from vorticity and divergence
!
call VorDiv2UV( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xyz_Vor = xyz_VorN, xyz_Div = xyz_DivN, & ! (in)
& xyz_U = xyz_UN, xyz_V = xyz_VN ) ! (out)
!---------------------------------------------
! 格子点上での非線形力学項の計算
! Calculate non-linear dynamical terms on grid points
!
call NonLinearOnGrid( dyn_as = dyn_sp_as % dyn_as, & ! (inout)
& xyz_U = xyz_UN, xyz_V = xyz_VN, & ! (in)
& xyz_Vor = xyz_VorN, xyz_Div = xyz_DivN, & ! (in)
& xyz_Temp = xyz_TempN, xyz_QVap = xyz_QVapN, & ! (in)
& xy_GradLonPi = xy_GradLonPiN, xy_GradLatPi = xy_GradLatPiN, & ! (in)
& xyz_UAdv = xyz_UAdvN, xyz_VAdv = xyz_VAdvN, & ! (out)
& xyz_DTempDt = xyz_DTempDtN, xyz_DQVapDt = xyz_DQVapDtN, & ! (out)
& xyz_KE = xyz_KEN, & ! (out)
& xyz_TempUAdv = xyz_TempUAdvN, xyz_TempVAdv = xyz_TempVAdvN, &
& xy_DPiDt = xy_DPiDtN, & ! (out)
& xyz_QVapUAdv = xyz_QVapUAdvN, xyz_QVapVAdv = xyz_QVapVAdvN ) ! (out)
!---------------------------------------------
! スペクトル時間変化項の計算
! Calculate spectral tendency terms
!
call Tendency( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xyz_UAdv = xyz_UAdvN, xyz_VAdv = xyz_VAdvN, xyz_KE = xyz_KEN, & ! (in)
& xyz_TempUAdv = xyz_TempUAdvN, & ! (in)
& xyz_TempVAdv = xyz_TempVAdvN, xyz_DTempDt = xyz_DTempDtN, & ! (in)
& xy_DPiDt = xy_DPiDtN, & ! (in)
& xyz_QVapUAdv = xyz_QVapUAdvN, xyz_QVapVAdv = xyz_QVapVAdvN, &
& xyz_DQVapDt = xyz_DQVapDtN, & ! (in)
& wz_DVorDt = wz_DVorDtN, wz_DDivDt = wz_DDivDtN, & ! (out)
& wz_DTempDt = wz_DTempDtN, w_DPiDt = w_DPiDtN, & ! (out)
& wz_DQVapDt = wz_DQVapDtN ) ! (out)
!---------------------------------------------
! 格子点値をスペクトル値へ ( $ t-\Delta t$ )
! Exchange grid values to spectral values ( $ t-\Delta t$ )
!
call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xyz_Vor = xyz_VorB, xyz_Div = xyz_DivB, xyz_Temp = xyz_TempB, & ! (in)
& xy_Ps = xy_PsB, xyz_QVap = xyz_QVapB, & ! (in)
& wz_Vor = wz_VorB, wz_Div = wz_DivB, wz_Temp = wz_TempB, & ! (out)
& w_Pi = w_PiB, wz_QVap = wz_QVapB ) ! (out)
!-----------------------------------------------------------------
! 時間積分
! Time integration
!-----------------------------------------------------------------
!---------------------------------------------
! セミインプリシット法による時間積分
! Time integration by semi-Implicit method
!
call TimeIntegration( dyn_as = dyn_sp_as % dyn_as, & ! (inout)
& wz_DVorDtN = wz_DVorDtN, wz_DDivDtN = wz_DDivDtN, & ! (in)
& wz_DTempDtN = wz_DTempDtN, wz_DQVapDtN = wz_DQVapDtN, & ! (in)
& w_DPiDtN = w_DPiDtN, & ! (in)
& wz_VorB = wz_VorB, wz_DivB = wz_DivB, & ! (in)
& wz_TempB = wz_TempB, wz_QVapB = wz_QVapB, & ! (in)
& w_PiB = w_PiB, & ! (in)
& wz_VorA = wz_VorA, wz_DivA = wz_DivA, & ! (out)
& wz_TempA = wz_TempA, wz_QVapA = wz_QVapA, & ! (out)
& w_PiA = w_PiA ) ! (out)
!---------------------------------------------
! スペクトル値を格子点値へ ( $ t+\Delta t$ )
! Exchange spectral values to grid values ( $ t+\Delta t$ )
!
call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& wz_Vor = wz_VorA, wz_Div = wz_DivA, wz_Temp = wz_TempA, & ! (in)
& w_Pi = w_PiA, wz_QVap = wz_QVapA, & ! (in)
& xyz_Vor = xyz_VorA, xyz_Div = xyz_DivA, xyz_Temp = xyz_TempA, & ! (out)
& xy_Ps = xy_PsA, xyz_QVap = xyz_QVapA, & ! (out)
& xyz_U = xyz_UA, xyz_V = xyz_VA ) ! (out)
!-----------------------------------------------------------------
! 拡散による補正
! Correction by diffusion
!-----------------------------------------------------------------
!---------------------------------------------
! 運動量水平拡散による渦度発散の時間変化
! Vorticity and divergence tendency by
! horizontal diffusion of momentum
!
call DiffusionVorDiv( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& wz_Vor = wz_VorA, wz_Div = wz_DivA, & ! (in)
& wz_VorDiff = wz_VorDiffA, wz_DivDiff = wz_DivDiffA ) ! (out)
!---------------------------------------------
! 運動量水平拡散による摩擦熱補正
! Frictional thermal correction by horizontal momentum diffusion
!
call DiffusionCorrectTemp( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& wz_VorDiff = wz_VorDiffA, wz_DivDiff = wz_DivDiffA, & ! (in)
& xyz_U = xyz_UA, xyz_V = xyz_VA, & ! (in)
& xyz_Temp = xyz_TempA ) ! (inout)
!-----------------------------------------------------------------
! 診断値の出力
! Output of diagnostic values
!-----------------------------------------------------------------
! ※ dyn_as の手続きを呼び出す. (SigmaDot などを出力予定).
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83Dynamics
subroutine DynSpectralAS83VorDiv2UV( dyn_sp_as, &
& xyz_Vor, xyz_Div, &
& xyz_U, xyz_V, &
& err )
!
! 与えられた渦度と発散から東西風速と南北風速を計算します.
!
! なお, 与えられた *dyn_sp_as* が Create によって初期設定
! されていない場合, プログラムはエラーを発生させます.
!
! Calculate zonal and meridional wind from given
! vorticity and divergence.
!
! If *dyn_sp_as* is not initialized by Create yet,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT, DP
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dyn_spectral, only: VorDiv2UV
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
real(DP), intent(in):: xyz_Vor (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 渦度 $ \zeta $
real(DP), intent(in):: xyz_Div (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 発散 $ D $
real(DP), intent(out):: xyz_U (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 東西風速 $ U $
real(DP), intent(out):: xyz_V (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 南北風速 $ V $
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83VorDiv2UV'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (.not. dyn_sp_as % initialized) then
stat = DCPAM_ENOTINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! 渦度発散から風速の計算
! Calculate wind velocity from vorticity and divergence
!-----------------------------------------------------------------
call VorDiv2UV( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xyz_Vor = xyz_Vor, xyz_Div = xyz_Div, & ! (in)
& xyz_U = xyz_U, xyz_V = xyz_V ) ! (out)
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83VorDiv2UV
subroutine DynSpectralAS83UV2VorDiv( dyn_sp_as, &
& xyz_Vor, xyz_Div, &
& xyz_U, xyz_V, &
& err )
!
! 与えられた東西風速と南北風速から渦度と発散を計算します.
!
! なお, 与えられた *dyn_sp_as* が Create によって初期設定
! されていない場合, プログラムはエラーを発生させます.
!
! Calculate vorticity and divergence
! from given zonal and meridional wind.
!
! If *dyn_sp_as* is not initialized by Create yet,
! error is occurred.
!
use dc_trace, only: BeginSub, EndSub
use dc_types, only: STRING, STDOUT, DP
use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
use dyn_spectral, only: UV2VorDiv
implicit none
type(DYNSPAS83), intent(inout):: dyn_sp_as
real(DP), intent(in):: xyz_U (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 東西風速 $ U $
real(DP), intent(in):: xyz_V (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 南北風速 $ V $
real(DP), intent(out):: xyz_Vor (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 渦度 $ \zeta $
real(DP), intent(out):: xyz_Div (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1)
! 発散 $ D $
logical, intent(out), optional:: err
! 例外処理用フラグ.
! デフォルトでは, この手続き内でエラーが
! 生じた場合, プログラムは強制終了します.
! 引数 *err* が与えられる場合,
! プログラムは強制終了せず, 代わりに
! *err* に .true. が代入されます.
!
! Exception handling flag.
! By default, when error occur in
! this procedure, the program aborts.
! If this *err* argument is given,
! .true. is substituted to *err* and
! the program does not abort.
integer:: stat
character(STRING):: cause_c
character(*), parameter:: subname = 'DynSpectralAS83UV2VorDiv'
continue
call BeginSub(subname)
stat = DC_NOERR
cause_c = ''
!-----------------------------------------------------------------
! 初期設定のチェック
! Check initialization
!-----------------------------------------------------------------
if (.not. dyn_sp_as % initialized) then
stat = DCPAM_ENOTINIT
cause_c = 'DYNSPAS83'
goto 999
end if
!-----------------------------------------------------------------
! 風速から渦度発散の計算
! Calculate vorticity and divergence from wind velocity
!-----------------------------------------------------------------
call UV2VorDiv( dyn_sp = dyn_sp_as % dyn_sp, & ! (inout)
& xyz_U = xyz_U, xyz_V = xyz_V, & ! (in)
& xyz_Vor = xyz_Vor, xyz_Div = xyz_Div ) ! (out)
!-----------------------------------------------------------------
! 終了処理, 例外処理
! Termination and Exception handling
!-----------------------------------------------------------------
999 continue
call StoreError(stat, subname, err, cause_c)
call EndSub(subname)
end subroutine DynSpectralAS83UV2VorDiv
!!$ subroutine DynSpectralAS83Sample( dyn_sp_as, err )
!!$ !--
!!$ ! DynSpectralAS83Sample の要約を記述してください.
!!$ !++
!!$ ! なお, 与えられた *dyn_sp_as* が Create によって初期設定
!!$ ! されていない場合, プログラムはエラーを発生させます.
!!$ !--
!!$ ! Describe brief of DynSpectralAS83Sample
!!$ !++
!!$ ! If *dyn_sp_as* is not initialized by Create yet,
!!$ ! error is occurred.
!!$ !
!!$ use dc_trace, only: BeginSub, EndSub
!!$ use dc_types, only: STRING, STDOUT
!!$ use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT
!!$ implicit none
!!$ type(DYNSPAS83), intent(inout):: dyn_sp_as
!!$ logical, intent(out), optional:: err
!!$ ! 例外処理用フラグ.
!!$ ! デフォルトでは, この手続き内でエラーが
!!$ ! 生じた場合, プログラムは強制終了します.
!!$ ! 引数 *err* が与えられる場合,
!!$ ! プログラムは強制終了せず, 代わりに
!!$ ! *err* に .true. が代入されます.
!!$ !
!!$ ! Exception handling flag.
!!$ ! By default, when error occur in
!!$ ! this procedure, the program aborts.
!!$ ! If this *err* argument is given,
!!$ ! .true. is substituted to *err* and
!!$ ! the program does not abort.
!!$ integer:: stat
!!$ character(STRING):: cause_c
!!$ character(*), parameter:: subname = 'DynSpectralAS83Sample'
!!$ continue
!!$ call BeginSub(subname)
!!$ stat = DC_NOERR
!!$ cause_c = ''
!!$
!!$ !-----------------------------------------------------------------
!!$ ! 初期設定のチェック
!!$ ! Check initialization
!!$ !-----------------------------------------------------------------
!!$ if (.not. dyn_sp_as % initialized) then
!!$ stat = DCPAM_ENOTINIT
!!$ cause_c = 'DYNSPAS83'
!!$ goto 999
!!$ end if
!!$
!!$
!!$
!!$ !-----------------------------------------------------------------
!!$ ! 終了処理, 例外処理
!!$ ! Termination and Exception handling
!!$ !-----------------------------------------------------------------
!!$999 continue
!!$ call StoreError(stat, subname, err, cause_c)
!!$ call EndSub(subname)
!!$ end subroutine DynSpectralAS83Sample
end module dyn_spectral_as83