Class x_base_module
In: setup/x_base_module.f90

1 次元 (x 方向) 等間隔交互格子 有限差分モデル用 基本モジュール

概要

x_base_module は, 1 次元 (x 方向) 等間隔交互格子を用いた有限差分 法に基づく数値モデルのための, 基本的な Fortran90 副プログラムおよび 関数を提供する.

このモジュールは x_module の下位モジュールである. 下請けモジュール として data_type モジュールを用いている.

備考

  • 例外処理がない
  • その場合のメッセージ表示がない
  • 初期化の値はマシンイプシロン値としておくべき

Methods

AvrX_p   AvrX_x   IntX_p   IntX_x   im   imax   imin   p_X   p_avr_x   p_dx   x_X   x_avr_p   x_axis_init   x_dx   xmargin  

Included Modules

dc_types

Public Instance methods

Function :
AvrX_p :real(DBKIND)
: 出力
p_Var(imin:imax) :real(DBKIND), intent(in)
: 入力

整数格子上の配列に対し x 方向の平均操作を行う

[Source]

    function AvrX_p(p_Var)
      ! 整数格子上の配列に対し x 方向の平均操作を行う

      real(DBKIND), intent(in) :: p_Var(imin:imax)  ! 入力
      real(DBKIND)             :: AvrX_p            ! 出力

      AvrX_p = IntX_p(p_Var)/sum(p_dx(1:im))

    end function AvrX_p
Function :
AvrX_x :real(DBKIND)
: 出力
x_Var(imin:imax) :real(DBKIND), intent(in)
: 入力

半整数格子上の配列に対し x 方向の平均操作を行う

[Source]

    function AvrX_x(x_Var)
      ! 半整数格子上の配列に対し x 方向の平均操作を行う

      real(DBKIND), intent(in) :: x_Var(imin:imax)  ! 入力
      real(DBKIND)             :: AvrX_x            ! 出力

      AvrX_x = IntX_x(x_Var)/sum(x_dx(1:im))

    end function AvrX_x
Function :
IntX_p :real(DBKIND)
: 出力
p_Var(imin:imax) :real(DBKIND), intent(in)
: 入力

整数格子上の配列に対し x 方向に重み付きの積分を行う

[Source]

    function IntX_p(p_Var)
      ! 整数格子上の配列に対し x 方向に重み付きの積分を行う

      real(DBKIND), intent(in) :: p_Var(imin:imax)  ! 入力
      real(DBKIND)             :: IntX_p            ! 出力

      ! 初期化
      IntX_p = 0.0d0

      ! 積分
      IntX_p = sum(p_Var(1:im)*p_dx(1:im))

    end function IntX_p
Function :
IntX_x :real(DBKIND)
: 出力
x_Var(imin:imax) :real(DBKIND), intent(in)
: 入力

半整数格子上の配列に対し x 方向に重み付きの積分を行う

[Source]

    function IntX_x(x_Var)
      ! 半整数格子上の配列に対し x 方向に重み付きの積分を行う

      real(DBKIND), intent(in) :: x_Var(imin:imax)  ! 入力
      real(DBKIND)             :: IntX_x            ! 出力
      real(DBKIND), allocatable:: work(:)

      ! 初期化
      IntX_x = 0.0d0
      allocate(work(imin:imax))

      ! 積分
      work = x_Var*x_dx

      IntX_x = sum(work(1:im))
      
      deallocate(work)

    end function IntX_x
im
Variable :
im = 10 :integer
: 格子点数
imax
Variable :
imax :integer
: 配列添字の上限
imin
Variable :
imin :integer
: 配列添字の下限
p_X
Variable :
p_X(:) :real(DBKIND),allocatable
: 整数格子点座標
Function :
p_avr_x(imin:imax) :real(DBKIND)
: 出力
x_Var(imin:imax) :real(DBKIND),intent(in)
: 入力

平均操作を行い半整数格子点の配列値を整数格子点上へ返す

[Source]

    function p_avr_x(x_Var)
      ! 平均操作を行い半整数格子点の配列値を整数格子点上へ返す
  
      real(DBKIND),intent(in) :: x_Var(imin:imax)  ! 入力
      real(DBKIND)            :: p_avr_x(imin:imax) ! 出力
      integer             :: ix                ! ループ添字

      ! 初期化
      ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき
      !
      p_avr_x = 0.0d0

      ! 平均操作
      ! * 平均がとれない imax 格子上の値は計算しない.
      ! * 不等間隔格子なので重みをつけて平均する.
      !   点 P, Q を a:b に内分する点の X の値は, a*X(Q) + b*X(P) 
      !   (a+b=1 を仮定) であることから, 
      !
      !     p_Var(ix) = [0.5*x_dx(ix)  /p_dx(ix)]*x_Var(ix+1) +    
      !                 [0.5*x_dx(ix+1)/p_dx(ix)]*x_Var(ix)     
      !   
      !   変形すると
      !
      !     p_Var(ix) = [x_dx(ix)*x_Var(ix+1) + x_dx(ix+1)*x_Var(ix)]
      !                 *0.5/p_dx(ix)
      ! 
      do ix = imin, imax-1
        p_avr_x(ix) = (x_dx(ix)*x_Var(ix+1) + x_dx(ix+1)*x_Var(ix))*0.5d0/p_dx(ix)
      end do

      ! imax 格子点上の値
      p_avr_x(imax) = p_avr_x(imax-1)

    end function p_avr_x
p_dx
Variable :
p_dx(:) :real(DBKIND),allocatable
: 整数格子点座標
x_X
Variable :
x_X(:) :real(DBKIND),allocatable
: 半整数格子点座標
Function :
x_avr_p(imin:imax) :real(DBKIND)
: 出力
p_Var(imin:imax) :real(DBKIND),intent(in)
: 入力

平均操作を行い整数格子点の配列値を半整数格子点上へ返す

[Source]

    function x_avr_p(p_Var)
      ! 平均操作を行い整数格子点の配列値を半整数格子点上へ返す
  
      real(DBKIND),intent(in) :: p_Var(imin:imax)  ! 入力
      real(DBKIND)            :: x_avr_p(imin:imax) ! 出力
      integer             :: ix                ! ループ添字

      ! 初期化
      ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき
      !
      x_avr_p = 0.0d0

      ! 平均操作
      ! * 平均がとれない imin 格子上の値は計算しない.
      !
      do ix = imin+1, imax
        x_avr_p(ix) = (p_Var(ix) + p_Var(ix-1))*0.5d0 
      end do

      ! imin 格子上の値
      x_avr_p(imin) = x_avr_p(imin+1)

    end function x_avr_p
Subroutine :
i :integer,intent(in)
: x 方向格子点数
xmg :integer,intent(in)
: x 方向糊代格子点数
xmin :real(DBKIND),intent(in)
: x 座標最小値
xmax :real(DBKIND),intent(in)
: x 座標最大値

配列の上下限の値と格子点座標値, 格子点間隔を設定する

[Source]

    subroutine x_axis_init(i, xmg, xmin, xmax)
      ! 配列の上下限の値と格子点座標値, 格子点間隔を設定する

      integer,intent(in)  :: i    ! x 方向格子点数
      integer,intent(in)  :: xmg  ! x 方向糊代格子点数
      real(DBKIND),intent(in) :: xmin ! x 座標最小値     
      real(DBKIND),intent(in) :: xmax ! x 座標最大値  
      integer             :: ix   ! do ループ添字
      real(DBKIND)            :: dx   ! 格子間隔

      im = i
      xmargin = xmg

      imin = 1  - xmargin
      imax = im + xmargin

      allocate(x_X(imin:imax))
      allocate(p_X(imin:imax))
      allocate(x_dx(imin:imax))
      allocate(p_dx(imin:imax))

      dx = (xmax - xmin)/im

      do ix = imin, imax
!        p_X(ix) = dx * ix  + myrank * dx * imax
!        x_X(ix) = dx * (ix - 0.5) + myrank * dx * imax
        p_X(ix) = dx * ix  
        x_X(ix) = dx * (ix - 0.5)
        x_dx(ix) = dx
        p_dx(ix) = dx
      end do

    end subroutine x_axis_init  
x_dx
Variable :
x_dx(:) :real(DBKIND),allocatable
: 半整数格子点座標
xmargin
Variable :
xmargin = 2 :integer
: 糊代格子点数