Class gridset
In: setup/gridset.f90

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

Methods

DelX   DelZ   DimXMax   DimXMin   DimZMax   DimZMin   FileNX   FileNZ   FileXMax   FileXMin   FileZMax   FileZMin   MarginX   MarginZ   NX   NZ   RegXMax   RegXMin   RegZMax   RegZMin   SpcNum   Xmax   Xmin   Zmax   Zmin   f_X   f_Z   gridset_init   s_X   s_Z  

Included Modules

DebugSet chemcalc

Public Instance methods

DelX()
Variable :
DelX :real(8)
: 格子間隔
DelZ()
Variable :
DelZ :real(8)
: 格子間隔
DimXMax()
Variable :
DimXMax :integer
: x 方向の配列の上限
DimXMin()
Variable :
DimXMin :integer
: x 方向の配列の下限
DimZMax()
Variable :
DimZMax :integer
: z 方向の配列の上限
DimZMin()
Variable :
DimZMin :integer
: z 方向の配列の下限
FileNX()
Variable :
FileNX :integer
: ファイル出力用
FileNZ()
Variable :
FileNZ :integer
: ファイル出力用
FileXMax()
Variable :
FileXMax :integer
: ファイル出力用
FileXMin()
Variable :
FileXMin :integer
: ファイル出力用
FileZMax()
Variable :
FileZMax :integer
: ファイル出力用
FileZMin()
Variable :
FileZMin :integer
: ファイル出力用
MarginX()
Constant :
MarginX = 5 :integer, parameter
: 境界のグリッド数
MarginZ()
Constant :
MarginZ = 5 :integer, parameter
: 境界のグリッド数
NX()
Variable :
NX :integer
: 格子点数
NZ()
Variable :
NZ :integer
: 格子点数
RegXMax()
Variable :
RegXMax :integer
: x 方向の物理領域の上限
RegXMin()
Variable :
RegXMin :integer
: x 方向の物理領域の下限
RegZMax()
Variable :
RegZMax :integer
: z 方向の物理領域の上限
RegZMin()
Variable :
RegZMin :integer
: z 方向の物理領域の下限
SpcNum()
Variable :
SpcNum :integer
: 化学種の数
Xmax()
Variable :
Xmax :real(8)
: x 座標の始点・終点
Xmin()
Variable :
Xmin :real(8)
: x 座標の始点・終点
Zmax()
Variable :
Zmax :real(8)
: z 座標の始点・終点
Zmin()
Variable :
Zmin :real(8)
: z 座標の始点・終点
f_X()
Variable :
f_X(:) :real(8), allocatable
: X 座標軸(ベクトル格子点)
f_Z()
Variable :
f_Z(:) :real(8), allocatable
: Z 座標軸(ベクトル格子点)
Subroutine :
cfgfile :character(*), intent(in)

NAMELIST から情報を得て, 格子点を計算する

This procedure input/output NAMELIST#gridset .

[Source]

  subroutine gridset_init(cfgfile)
    !
    !NAMELIST から情報を得て, 格子点を計算する
    !

    !モジュール読み込み
    use DebugSet, only: DebugOn
    use chemcalc, only: ChemCalcDim_init !配列の大きさの初期化
    
    !暗黙の型宣言禁止
    implicit none

    !入力変数
    character(*), intent(in) :: cfgfile

    !内部変数
    integer            :: i, k
    integer, parameter :: kind = 8      !精度を表す

    !-----------------------------------------------------------------
    ! NAMELIST から情報を取得
    !-----------------------------------------------------------------
    NAMELIST /gridset/ NX, NZ, Xmin, Xmax, Zmin, Zmax, SpcNum

    open (10, FILE=cfgfile)
    read(10, NML=gridset)
    close(10)

    !-----------------------------------------------------------------
    ! 格子点間隔計算
    !-----------------------------------------------------------------
    DelX = (Xmax - Xmin) / real(NX, kind)
    DelZ = (Zmax - Zmin) / real(NZ, kind)

    !-----------------------------------------------------------------
    ! 物理的に意味のある領域の上限・下限を決める
    !-----------------------------------------------------------------;
    RegXMin = 0
    RegXMax = NX
    RegZMin = 0
    RegZMax = NZ

    !-----------------------------------------------------------------
    ! 配列の上限・下限を決める
    !-----------------------------------------------------------------
    DimXMin = RegXMin - MarginX
    DimXMax = RegXMax + MarginX
    DimZMin = RegZMin - MarginZ
    DimZMax = RegZMax + MarginZ

    !-----------------------------------------------------------------
    ! グリッドの設定
    !-----------------------------------------------------------------
    allocate( f_X(DimXMin:DimXMax), f_Z(DimZMin:DimZMax) )
    do i = DimXMin, DimXMax
      f_X(i) = Xmin + DelX * ( real(i, kind) - RegXMin )
    end do
    do k = DimZMin, DimZMax
      f_Z(k) = Zmin + DelZ * ( real(k, kind) - RegZMin )
    end do
    
    !-----------------------------------------------------------------
    ! 半格子ずれたグリッドの設定. 等間隔を前提にしている.
    !-----------------------------------------------------------------
    allocate( s_X(DimXMin: DimXMax), s_Z(DimZMin: DimZMax) )
    s_X(DimXMin+1:DimXMax) = ( f_X(DimXMin+1:DimXMax) + f_X(DimXMin:DimXMax-1) ) * 5.0d-1
    s_X(DimXMin) = s_X(DimXMin+1) - DelX
    
    s_Z(DimZMin+1:DimZMax) = ( f_Z(DimZMin+1:DimZMax) + f_Z(DimZMin:DimZMax-1) ) * 5.0d-1
    s_Z(DimZMin) = s_Z(DimZMin+1) - DelZ

    !-----------------------------------------------------------------
    ! デバッグモードか否かで, ファイル出力に用いる変数の大きさを変える
    !-----------------------------------------------------------------
    if (DebugOn) then
      FileNX   = size(f_X, 1)
      FileNZ   = size(f_Z, 1)
      FileXMin = DimXMin
      FileXMax = DimXMax
      FileZMin = DimZMin
      FileZMax = DimZMax
    else
      FileNX   = NX
      FileNZ   = NZ
      FileXMin = RegXMin + 1
      FileXMax = RegXMax
      FileZMin = RegZMin + 1
      FileZMax = RegZMax
    end if

    !-----------------------------------------------------------------    
    ! 確認
    !-----------------------------------------------------------------
    write(*,*) "gridset_init, XMin   ", XMin
    write(*,*) "gridset_init: XMax   ", XMax
    write(*,*) "gridset_init: ZMin   ", ZMin
    write(*,*) "gridset_init: ZMax   ", ZMax
    write(*,*) "gridset_init: DelX   ", DelX
    write(*,*) "gridset_init: DelZ   ", DelZ
    write(*,*) "gridset_init: NX     ", NX
    write(*,*) "gridset_init: NZ     ", NZ
    write(*,*) "gridset_init: SpcNum ", SpcNum
    write(*,*) "gridset_init: MarginX", MarginX
    write(*,*) "gridset_init: MarginZ", MarginZ
    write(*,*) "gridset_init: DimXMin", DimXMin
    write(*,*) "gridset_init: DimXMax", DimXMax
    write(*,*) "gridset_init: DimZMin", DimZMin
    write(*,*) "gridset_init: DimZMax", DimZMax
  end subroutine gridset_init
s_X()
Variable :
s_X(:) :real(8), allocatable
: X 座標軸(スカラー格子点)
s_Z()
Variable :
s_Z(:) :real(8), allocatable
: Z 座標軸(スカラー格子点)

[Validate]