Class gridset
In: setup/gridset.F90

格子点数・最大波数設定

Number of grid points and maximum truncated wavenumber settings

Note that Japanese and English are described in parallel.

格子点数の設定および保管を行います. スペクトル法を用いることを前提にしており, 最大波数の設定と保管も行います.

Number of grid points is set and stored. Maximum truncated wavenumber is set and stored too, because spectral method is expected to be used.

Variables List

nmax :最大全波数
lmax :スペクトルデータの配列サイズ
imax :経度格子点数
jmax_global :緯度格子点数 (全球)
jmax :緯度格子点数
kmax :鉛直層数
kslmax :地下鉛直層数
———— :————
nmax :Maximum truncated wavenumber
lmax :Size of array for spectral data
imax :Number of grid points in longitude
jmax_global :Number of grid points in latitude on whole globe
jmax :Number of grid points in latitude
kmax :Number of vertical level
kslmax :Number of subsurface vertical level

Procedures List

GridsetInit :格子点数と最大波数の設定
———— :————
GridsetInit :Settings of number of grid points and maximum truncated wavenumber

NAMELIST

NAMELIST#gridset_nml

Methods

Included Modules

dc_types namelist_util dc_iounit dc_message mpi_wrapper

Public Instance methods

Subroutine :
jc :integer, intent(in )

それぞれのプロセスが担当する緯度格子点数を確認.

Check number of latitudinal grids on each process.

[Source]

  subroutine GridsetCheckNumberOfLatGrid( jc )
    !
    ! それぞれのプロセスが担当する緯度格子点数を確認.
    !
    ! Check number of latitudinal grids on each process.
    !

    ! モジュール引用 ; USE statements
    !

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! MPI
    !
    use mpi_wrapper, only: nprocs, myrank, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait


    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer, intent(in ) :: jc

    ! Local variable
    !
    integer , allocatable :: a_jmax_tmp(:)
    real(DP)              :: a_SendBuf (1)
    real(DP), allocatable :: aa_RecvBuf(:,:)
    integer , allocatable :: a_iReqSend(:)
    integer , allocatable :: a_iReqRecv(:)
    integer               :: n

    ! 実行文 ; Executable statement
    !


    if ( jmax /= jc ) then
      n = myrank
      call MessageNotify( 'E', module_name, 'Saved jmax in myproc=%d is %d, but is not the same as one calculated by spml, %d.', i = (/myrank, jmax, jc/) )
    end if



    allocate( a_jmax_tmp(0:nprocs-1) )
    allocate( aa_RecvBuf(1,0:nprocs-1) )
    allocate( a_iReqSend(0:nprocs-1) )
    allocate( a_iReqRecv(0:nprocs-1) )

    a_SendBuf(1) = jc

    ! jmax is sent to all nodes, and received by all nodes.

    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      call MPIWrapperISend( n, 1, a_SendBuf      , a_iReqSend(n) )
      call MPIWrapperIRecv( n, 1, aa_RecvBuf(:,n), a_iReqRecv(n) )
    end do

    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      call MPIWrapperWait( a_iReqSend(n) )
      call MPIWrapperWait( a_iReqRecv(n) )
    end do


    do n = 0, nprocs-1
      if ( n == myrank ) then
        a_jmax_tmp(n) = jc
      else
        a_jmax_tmp(n) = aa_RecvBuf(1,n)
      end if
    end do

    ! Check number of latitudinal grid in each process.
    do n = 0, nprocs-1
      if ( a_jmax(n) /= a_jmax_tmp(n) ) then
        call MessageNotify( 'E', module_name, 'jmax in proc=%d is %d, but is not the same as saved one, %d.', i = (/n, a_jmax_tmp(n), a_jmax(n)/) )
      end if
    end do


    deallocate( aa_RecvBuf )
    deallocate( a_iReqSend )
    deallocate( a_iReqRecv )
    deallocate( a_jmax_tmp )



  end subroutine GridsetCheckNumberOfLatGrid
Subroutine :

gridset モジュールの初期化を行います. NAMELIST#gridset_nml の読み込みはこの手続きで行われます.

"gridset" module is initialized. NAMELIST#gridset_nml is loaded in this procedure.

This procedure input/output NAMELIST#gridset_nml .

[Source]

  subroutine GridsetInit
    !
    ! gridset モジュールの初期化を行います. 
    ! NAMELIST#gridset_nml の読み込みはこの手続きで行われます. 
    !
    ! "gridset" module is initialized. 
    ! NAMELIST#gridset_nml is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 宣言文 ; Declaration statements
    !
    implicit none
    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /gridset_nml/ nmax, imax, jmax, kmax, kslmax
          !
          ! デフォルト値については初期化手続 "gridset#GridsetInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "gridset#GridsetInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( gridset_inited ) return
    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
#ifdef AXISYMMETRY
    imax        = 1
#elif AXISYMMETRY_SJPACK
    imax        = 1
#else
    imax        = 32
#endif
    jmax        = imax / 2
    kmax        = 5
    nmax        = ( imax - 1 ) / 3

    kslmax      = 0

    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = gridset_nml, iostat = iostat_nml ) ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
      if ( iostat_nml == 0 ) write( STDOUT, nml = gridset_nml )
    end if

    ! Check of array size
    !
#ifdef AXISYMMETRY
    if ( imax /= 1 ) then
      call MessageNotify( 'E', module_name, 'Number of grid points in longitudinal direction has to be 1.' )
    end if
#elif AXISYMMETRY_SJPACK
    if ( imax /= 1 ) then
      call MessageNotify( 'E', module_name, 'Number of grid points in longitudinal direction has to be 1.' )
    end if
#endif

    if ( imax /= 1 ) then
      if ( jmax /= imax / 2 ) then
        call MessageNotify( 'E', module_name, 'Number of grid points in latitudinal direction has to be equal to ' // 'number of grid points in longitudinal direction divided by 2.' )
      end if
    end if

    ! Set number of grid points in latitudinal direction
    !
    jmax_global = jmax
    !
    ! Calculate and save jmax in all processes
    !   If jmax is 1, the model runs as 1D model. In that case, jmax need not 
    !   be calculated. 
    ! Even if jmax has a value, it is overwritten. 
    !
    call GridsetCalcjmax


    ! Set size of array for spectral data
    !
#ifdef AXISYMMETRY
    lmax =  nmax+1
#elif AXISYMMETRY_SJPACK
    lmax =  nmax+1
#else
    lmax = (nmax+1)**2
#endif


    ! 格子点数のチェック
    ! Check number of grid points
    !
    if ( nmax < 1 .or. imax < 1 .or. jmax_global < 1 .or. jmax < 1 .or. kmax < 1 ) then
      call MessageNotify( 'E', module_name, 'number of grid points and maximum truncated wavenumber must be more than 0. ' // 'nmax=%d, imax=%d, jmax_global=%d, jmax=%d, kmax=%d' , i = (/ nmax, imax, jmax_global, jmax, kmax /) )
    end if

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )
    call MessageNotify( 'M', module_name, '    nmax        = %d', i = (/   nmax  /) )
    call MessageNotify( 'M', module_name, '    imax        = %d', i = (/   imax  /) )
    call MessageNotify( 'M', module_name, '    jmax_global = %d', i = (/   jmax_global /) )
    call MessageNotify( 'M', module_name, '    jmax        = %d', i = (/   jmax  /) )
    call MessageNotify( 'M', module_name, '    kmax        = %d', i = (/   kmax  /) )
    call MessageNotify( 'M', module_name, '  kslmax        = %d', i = (/ kslmax  /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    gridset_inited = .true.
  end subroutine GridsetInit
a_jmax
Variable :
a_jmax(:) :integer, save, allocatable, public
: 各プロセスの緯度格子点数. Number of grid points in latitude of each process
gridset_inited
Variable :
gridset_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag
imax
Variable :
imax :integer, save, public
: 経度格子点数. Number of grid points in longitude
jmax
Variable :
jmax :integer, save, public
: 緯度格子点数. Number of grid points in latitude
jmax_global
Variable :
jmax_global :integer, save, public
: 緯度格子点数 (全球). Number of grid points in latitude on whole globe
jmax_max
Variable :
jmax_max :integer, save , public
: 各プロセスの緯度格子点数の最大値. Maximum number of grid points in latitude of all processes
kmax
Variable :
kmax :integer, save, public
: 鉛直層数. Number of vertical level
kslmax
Variable :
kslmax :integer, save, public
: 地下の鉛直層数. Number of subsurface vertical level
lmax
Variable :
lmax :integer, save, public
: スペクトルデータの配列サイズ Size of array for spectral data
nmax
Variable :
nmax :integer, save, public
: 最大全波数. Maximum truncated wavenumber

Private Instance methods

Subroutine :

jmax を計算.

Calculate jmax

[Source]

  subroutine GridsetCalcjmax
    !
    ! jmax を計算.
    !
    ! Calculate jmax
    !

    ! モジュール引用 ; USE statements
    !

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! MPI
    !
    use mpi_wrapper, only: nprocs, myrank

    ! 宣言文 ; Declaration statements
    !
    implicit none


    ! Local variables
    !
    integer            :: jh, jph, js, je, jch, jc
    integer            :: n


    allocate( a_jmax(0:nprocs-1) )

    if ( ( imax == 1 ) .and. ( jmax == 1 ) ) then

      jmax_max = 1
      a_jmax   = 1

    else

      ! Following calculation is a copy of that in a subroutine SNMINI in snmini.f 
      ! included in ISPACK.

      jh = jmax_global / 2
      jph = ( jh - 1 ) / nprocs + 1

      jmax_max = jph * 2

      do n = 0, nprocs-1

        js = jph * n + 1
        je = min( jph * ( n + 1 ), jh )
        if ( je >= js ) THEN
          jch = je - js + 1
          jc  = jch * 2
        else
          jc = 0
          js = 1
          je = 1
        end if

        a_jmax(n) = jc

      end do

      if ( jmax /= a_jmax( myrank ) ) then
        call MessageNotify( 'M', module_name, 'jmax is overwritten, and is %d.', i = (/ a_jmax( myrank ) /) )
      end if
      jmax = a_jmax( myrank )

    end if


  end subroutine GridsetCalcjmax
Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_util_inited

    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

  end subroutine InitCheck
module_name
Constant :
module_name = ‘gridset :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20140204-3 $’ // ’$Id: gridset.F90,v 1.4 2012-02-01 12:03:32 yot Exp $’ :character(*), parameter
: モジュールのバージョン Module version