Class gridset
In: ../src/setup/gridset.f90

格子設定用モジュール

引数に与えられた NAMELIST ファイルから, 格子点情報を取得し, 保管するための変数参照型モジュール

Procedures List

gridset_init :初期化ルーチン
gridset_check :設定内容の整合性をチェックするためのルーチン

Methods

FlagCalc3D   NCMAX   NX   NY   NZ   Xmg   Ymg   Zmg   gridset_init   imax   imin   jmax   jmin   kmax   kmin   xsub   ysub  

Included Modules

dc_iounit dc_message mpi_wrapper namelist_util

Public Instance methods

FlagCalc3D
Variable :
FlagCalc3D = .true. :logical, public, save
NCMAX
Variable :
NCMAX = 1 :integer, public, save
: 組成配列要素数
NX
Variable :
NX = 10 :integer, public, save
: x 方向格子点数
NY
Variable :
NY = 10 :integer, public, save
: y 方向格子点数
NZ
Variable :
NZ = 10 :integer, public, save
: z 方向格子点数
Xmg
Variable :
Xmg = 2 :integer, public, save
: x 方向糊代格子点数
Ymg
Variable :
Ymg = 2 :integer, public, save
: y 方向糊代格子点数
Zmg
Variable :
Zmg = 2 :integer, public, save
: z 方向糊代格子点数
Subroutine :

初期化ルーチン

設定ファイルから情報を読み込み格子点数を計算する

This procedure input/output NAMELIST#gridset_nml .

[Source]

  subroutine gridset_init
    !
    != 初期化ルーチン
    !
    ! 設定ファイルから情報を読み込み格子点数を計算する
    !

    !暗黙の型宣言禁止
    implicit none

    !内部変数
    integer            :: unit                !設定ファイル用装置番号

    !-----------------------------------------------------------------
    ! 設定ファイルから情報を読み込み
    !
    NAMELIST /gridset_nml/ xdim, ydim, zdim, NCMAX, Xmg, Ymg, Zmg, xsub, ysub

    call FileOpen(unit, file=namelist_filename, mode='r')
    read(unit, NML=gridset_nml)
    close(unit)

    !-----------------------------------------------------------------    
    ! NX, NY, NZ を決める
    !
    nx = xdim / xsub 
    imin = 1  - xmg
    imax = nx + xmg
   
    if (ydim == 1) then
      ! 二次元
      ymg = 0
      ny = 1
      FlagCalc3D = .false. 
    else
      ! 三次元
      ny = ydim / ysub
      FlagCalc3D = .true. 
    end if
    jmin = 1  - ymg
    jmax = ny + ymg

    nz = zdim
    kmin = 1  - zmg
    kmax = nz + zmg    

    !-----------------------------------------------------------------    
    !"myrank == 0" に該当する計算ノードが, 読み込んだ情報を出力
    !
    if (myrank == 0) then 
      call MessageNotify( "M", "gridset_init", "xsub = %d",  i=(/xsub/) )
      call MessageNotify( "M", "gridset_init", "ysub = %d",  i=(/ysub/) )
      call MessageNotify( "M", "gridset_init", "xdim = %d",   i=(/xdim/) )
      call MessageNotify( "M", "gridset_init", "ydim = %d",   i=(/ydim/) )
      call MessageNotify( "M", "gridset_init", "zdim = %d",   i=(/zdim/) )
      call MessageNotify( "M", "gridset_init", "[1 node] NX = %d",   i=(/NX/) )
      call MessageNotify( "M", "gridset_init", "[1 node] NY = %d",   i=(/NY/) )
      call MessageNotify( "M", "gridset_init", "[1 node] NZ = %d",   i=(/NZ/) )
      call MessageNotify( "M", "gridset_init", "[1 node] NCMAX = %d",   i=(/NCMAX/) )
      call MessageNotify( "M", "gridset_init", "[1 node] xmg  = %d", i=(/Xmg/) )
      call MessageNotify( "M", "gridset_init", "[1 node] ymg  = %d", i=(/Ymg/) )
      call MessageNotify( "M", "gridset_init", "[1 node] zmg  = %d", i=(/Zmg/) )
      call MessageNotify( "M", "gridset_init", "[1 node] imin = %d", i=(/imin/) )
      call MessageNotify( "M", "gridset_init", "[1 node] imax = %d", i=(/imax/) )
      call MessageNotify( "M", "gridset_init", "[1 node] jmin = %d", i=(/jmin/) )
      call MessageNotify( "M", "gridset_init", "[1 node] jmax = %d", i=(/jmax/) )
      call MessageNotify( "M", "gridset_init", "[1 node] kmin = %d", i=(/kmin/) )
      call MessageNotify( "M", "gridset_init", "[1 node] kmax = %d", i=(/kmax/) )
    end if

    !-----------------------------------------------------------------
    ! 値のチェック
    !    
    call gridset_check

  end subroutine gridset_init
imax
Variable :
imax :integer, public, save
: x 方向の配列の上限
imin
Variable :
imin :integer, public, save
: x 方向の配列の下限
jmax
Variable :
jmax :integer, public, save
: y 方向の配列の下限
jmin
Variable :
jmin :integer, public, save
: y 方向の配列の下限
kmax
Variable :
kmax :integer, public, save
: z 方向の配列の下限
kmin
Variable :
kmin :integer, public, save
: z 方向の配列の下限
xsub
Variable :
xsub = 1 :integer, public, save
: Number of MPI processes (X direction)
ysub
Variable :
ysub = 1 :integer, public, save
: Number of MPI processes (Y direction)