| Class | et_module |
| In: |
src/et_module.f90
|
2 次元水路領域問題, Fourier 展開 + Chebyshev 展開法 spml/et_module モジュールは 2 次元水路領域での流体運動を スペクトル法により数値計算を実行するための Fortran90 関数を提供する. 周期的な境界条件を扱うための X 方向へのフーリエ変換と境界壁を扱うための Y 方向のチェビシェフ変換を用いる場合のスペクトル計算のためのさまざまな 関数を提供する. 内部で ae_module, at_module を用いている. 最下部ではフーリエ変換およびチェビシェフ変換のエンジンとして ISPACK/FTPACK の Fortran77 サブルーチンを用いている.
関数・変数の名前と型について
○命名法
* 関数名の先頭 (et_, yx_, x_, y_) は, 返す値の形を示している.
et_ : 2次元スペクトルデータ
yx_ : 2 次元格子点データ
x_ : X 方向 1 次元格子点データ
y_ : Y 方向 1 次元格子点データ
* 関数名の間の文字列(Dx, Dy, Lapla, LaplaInv, Jacobian)は,
その関数の作用を表している.
* 関数名の最後 (_et_et,_et,_yx, _x, _y) は, 入力変数のスペクトルデータ
および格子点データであることを示している.
_et : 2次元スペクトルデータ
_et_et : 2 つの2次元スペクトルデータ
_yx : 2 次元格子点データ
_x : X 方向 1 次元格子点データ
_y : Y 方向 1 次元格子点データ
○各データの種類の説明
* yx : 2 次元格子点データ.
変数の種類と次元は real(8), dimension(0:jm,0:im-1).
im, jm はそれぞれ X, Y 座標の格子点数であり, サブルーチン
et_initial にてあらかじめ設定しておく.
第 1 次元が Y 座標の格子点位置番号, 第 2 次元が X 座標の
格子点位置番号である (X, Y の順ではない)ことに注意.
* et : 2 次元スペクトルデータ.
変数の種類と次元は real(8), dimension(-km:km,0:lm).
km, lm はそれぞれ X, Y 方向の最大波数であり, サブルーチン
et_initial にてあらかじめ設定しておく.
スペクトルデータの格納のされ方については...
* x, y : X, Y 方向 1 次元格子点データ.
変数の種類と次元はそれぞれ real(8), dimension(0:im-1)
および real(8), dimension(0:jm).
* e, t : 1 次元スペクトルデータ.
変数の種類と次元は real(8), dimension(-km:km)
および real(8), dimension(-lm:lm).
* ax, ay : 1 次元格子点データの並んだ 2 次元配列.
変数の種類と次元は real(8), dimension(:,0:im-1)
および real(8), dimension(:,0:jm).
* ae, at : 1 次元スペクトルデータの並んだ 2 次元配列.
変数の種類と次元は real(8), dimension(:,-km:km)
および real(8), dimension(:,0:lm).
* et_ で始まる関数が返す値はスペクトルデータに同じ.
* yx_ で始まる関数が返す値は 2 次元格子点データに同じ.
* x_, y_ で始まる関数が返す値は 1 次元格子点データに同じ.
* スペクトルデータに対する微分等の作用とは, 対応する格子点データに
微分などを作用させたデータをスペクトル変換したものことである.
| Function : | |||
| AvrX_x : | real(8)
| ||
| x : | real(8), dimension(0:im-1)
|
1 次元(X)格子点データの X 方向平均
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.
function AvrX_x(x)
!
! 1 次元(X)格子点データの X 方向平均
!
! 実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し,
! x_X_Weight の総和で割ることで平均している.
!
real(8), dimension(0:im-1) :: x !(in) 1 次元格子点データ
real(8) :: AvrX_x !(out) 平均値
AvrX_x = IntX_x(x)/sum(x_X_weight)
end function AvrX_x
| Function : | |||
| AvrYX_yx : | real(8)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
2 次元格子点データの全領域平均
実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算し, x_X_Weight*y_Y_Weight の総和で割ることで平均している.
function AvrYX_yx(yx)
!
! 2 次元格子点データの全領域平均
!
! 実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた
! 総和を計算し, x_X_Weight*y_Y_Weight の総和で割ることで平均している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8) :: AvrYX_yx
!(out) 平均値
AvrYX_yx = IntYX_yx(yx)/(sum(x_X_weight)*sum(y_Y_weight))
end function AvrYX_yx
| Function : | |||
| AvrY_y : | real(8)
| ||
| y : | real(8), dimension(0:jm)
|
1 次元(Y)格子点データの Y 方向平均
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.
function AvrY_y(y)
!
! 1 次元(Y)格子点データの Y 方向平均
!
! 実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し,
! y_Y_Weight の総和で割ることで平均している.
!
real(8), dimension(0:jm) :: y !(in) 1 次元格子点データ
real(8) :: AvrY_y !(out) 平均値
AvrY_y = IntY_y(y)/sum(y_Y_weight)
end function AvrY_y
| Function : | |||
| IntX_x : | real(8)
| ||
| x : | real(8), dimension(0:im-1)
|
X 方向積分
1 次元(X)格子点データの X 方向積分
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
function IntX_x(x) ! X 方向積分
!
! 1 次元(X)格子点データの X 方向積分
!
! 実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
!
real(8), dimension(0:im-1) :: x !(in) 1 次元格子点データ
real(8) :: IntX_x !(out) 積分値
IntX_x = sum(x*x_X_Weight)
end function IntX_x
| Function : | |||
| IntYX_yx : | real(8)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
全領域積分
2 次元格子点データの全領域積分および平均.
実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算している.
function IntYX_yx(yx) ! 全領域積分
!
! 2 次元格子点データの全領域積分および平均.
!
! 実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた
! 総和を計算している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8) :: IntYX_yx
!(out) 積分値
integer :: i, j
IntYX_yx = 0.0d0
do i=0,im-1
do j=0,jm
IntYX_yx = IntYX_yx + yx(j,i) * y_Y_Weight(j) * x_X_Weight(i)
enddo
enddo
end function IntYX_yx
| Function : | |||
| IntY_y : | real(8)
| ||
| y : | real(8), dimension(0:jm)
|
Y 方向積分
1 次元(Y)格子点データの Y 方向積分
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
function IntY_y(y) ! Y 方向積分
!
! 1 次元(Y)格子点データの Y 方向積分
!
! 実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
!
real(8), dimension(0:jm) :: y !(in) 1 次元格子点データ
real(8) :: IntY_y !(out) 積分値
IntY_y = sum(y*y_Y_Weight)
end function IntY_y
| Function : | |||
| ae : | real(8), dimension(:,-km:), intent(in)
|
入力スペクトルデータに X 微分を作用する(2 次元データ).
スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.
The entity is ae_module#ae_Dx_ae
| Function : | |||
| ae_ag : | real(8), dimension(size(ag,1),-km:km)
| ||
| ag : | real(8), dimension(:,:), intent(in)
|
格子点データからスペクトルデータへ変換する(2 次元データ用)
The entity is ae_module#ae_ag
| Subroutine : | |||
| at_data : | real(8), dimension(:,0:),intent(inout)
| ||
| values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界での値を与える.
The entity is at_module#at_Boundaries_DD
| Subroutine : | |||
| t_data : | real(8), dimension(0:km),intent(inout)
| ||
| values : | real(8), dimension(2), intent(in), optional
|
Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界での値を与える.
The entity is at_module#at_Boundaries_DD
| Subroutine : | |||
| at_data : | real(8), dimension(:,0:),intent(inout)
| ||
| values : | real(8), dimension(:,:), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 2 次元配列用) i=0 で値, i=im で勾配の値を与える.
The entity is at_module#at_Boundaries_DN
| Subroutine : | |||
| t_data : | real(8), dimension(0:km),intent(inout)
| ||
| values : | real(8), dimension(2), intent(in), optional
|
Dirichlet/Neumann 型境界条件の適用(タウ法, 1 次元配列用) i=0 で値, i=im で勾配の値を与える.
The entity is at_module#at_Boundaries_DN
| Subroutine : | |||
| at_data : | real(8), dimension(:,0:),intent(inout)
| ||
| values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) i=0 で勾配の値, i=im で値を与える.
The entity is at_module#at_Boundaries_ND
| Subroutine : | |||
| t_data : | real(8), dimension(0:km),intent(inout)
| ||
| values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) i=0 で勾配の値, i=im で値を与える.
The entity is at_module#at_Boundaries_ND
| Subroutine : | |||
| at_data : | real(8), dimension(:,0:),intent(inout)
| ||
| values : | real(8), dimension(:,:), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 2 次元配列用) 両境界で勾配の値を与える.
The entity is at_module#at_Boundaries_NN
| Subroutine : | |||
| t_data : | real(8), dimension(0:km),intent(inout)
| ||
| values : | real(8), dimension(2), intent(in), optional
|
Neumann/Dirichlet 型境界条件の適用(タウ法, 1 次元配列用) 両境界で勾配の値を与える.
このサブルーチンを直接使うことを勧めない. 共通インターフェース at_Boundaries_NN を用いること.
The entity is at_module#at_Boundaries_NN
| Function : | |||
| at_Dx_at : | real(8), dimension(size(at_data,1),0:size(at_data,2)-1)
| ||
| at_data : | real(8), dimension(:,0:), intent(in)
|
入力チェビシェフデータに X 微分を作用する(2 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
The entity is at_module#at_Dx_at
| Function : | |||
| at_ag : | double precision, dimension(size(ag_data,1),0:km)
| ||
| ag_data : | double precision, dimension(:,:), intent(in)
|
格子データからチェビシェフデータへ変換する(2 次元配列用).
The entity is at_module#at_ag
| Function : | |||
| ag_ae : | real(8), dimension(size(ae,1),0:im-1)
| ||
| ae : | real(8), dimension(:,-km:), intent(in)
|
スペクトルデータから格子点データへ変換する(2 次元データ用)
The entity is ae_module#ag_ae
| Function : | |||
| ag_at : | double precision, dimension(size(at_data,1),0:im)
| ||
| at_data : | double precision, dimension(:,:), intent(in)
|
チェビシェフデータから格子データへ変換する(2 次元配列用).
The entity is at_module#ag_at
| Function : | |||
| e_Dx_e : | real(8), dimension(-km:km)
| ||
| e : | real(8), dimension(-km:km), intent(in)
|
入力スペクトルデータに X 微分を作用する(1 次元データ).
スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.
The entity is ae_module#e_Dx_e
| Function : | |||
| e_g : | real(8), dimension(-km:km)
| ||
| g : | real(8), dimension(0:im-1), intent(in)
|
格子 -> スペクトル
格子点データからスペクトルデータへ変換する(1 次元データ用)
The entity is ae_module#e_g
| Subroutine : | |||
| et : | real(8), dimension(-km:km,0:lm),intent(inout)
| ||
| values : | real(8), dimension(-km:km,2), intent(in), optional
| ||
| cond : | character(len=2), intent(in), optional
|
ディリクレ, ノイマン条件の適用. チェビシェフ空間での計算
実際には中で呼ばれている at_module のサブルーチン at_Boundaries_DD, at_Boundaries_DN, at_Boundaries_ND, at_Boundaries_NN を用いている. これらを直接呼ぶことも出来る.
subroutine et_Boundaries(et,values,cond)
!
! ディリクレ, ノイマン条件の適用. チェビシェフ空間での計算
!
! 実際には中で呼ばれている at_module のサブルーチン at_Boundaries_DD,
! at_Boundaries_DN, at_Boundaries_ND, at_Boundaries_NN を用いている.
! これらを直接呼ぶことも出来る.
!
real(8), dimension(-km:km,0:lm),intent(inout) :: et
! 境界条件を適用するデータ. 修正された値を返す.
real(8), dimension(-km:km,2), intent(in), optional :: values
! 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
! 省略時は値/勾配 0 となる.
character(len=2), intent(in), optional :: cond
! 境界条件. 省略時は 'DD'
! DD : 両端ディリクレ
! DN,ND : ディリクレ/ノイマン条件
! NN : 両端ノイマン
if (.not. present(cond)) then
if (present(values)) then
call at_Boundaries_DD(et,values)
else
call at_Boundaries_DD(et)
endif
return
endif
select case(cond)
case ('NN')
if (present(values)) then
call at_Boundaries_NN(et,values)
else
call at_Boundaries_NN(et)
endif
case ('DN')
if (present(values)) then
call at_Boundaries_DN(et,values)
else
call at_Boundaries_DN(et)
endif
case ('ND')
if (present(values)) then
call at_Boundaries_ND(et,values)
else
call at_Boundaries_ND(et)
endif
case ('DD')
if (present(values)) then
call at_Boundaries_DD(et,values)
else
call at_Boundaries_DD(et)
endif
case default
call MessageNotify('E','et_Boundaries','B.C. not supported')
end select
end subroutine et_Boundaries
| Function : | |
| et_Dx_et : | real(8), dimension(-km:km,0:lm) |
| et : | real(8), dimension(-km:km,0:lm), intent(in) |
入力スペクトルデータに X 微分(∂x)を作用する.
スペクトルデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのスペクトル変換のことである.
実際にはスペクトルデータに X 方向波数 k をかけて sin(kx) <-> cos(kx) 成分に入れ換える計算を行っている.
function et_Dx_et(et)
!
! 入力スペクトルデータに X 微分(∂x)を作用する.
!
! スペクトルデータの X 微分とは, 対応する格子点データに X 微分を
! 作用させたデータのスペクトル変換のことである.
!
! 実際にはスペクトルデータに X 方向波数 k をかけて
! sin(kx) <-> cos(kx) 成分に入れ換える計算を行っている.
!
real(8), dimension(-km:km,0:lm) :: et_Dx_et
real(8), dimension(-km:km,0:lm), intent(in) :: et
integer k
do k=-km,km
et_Dx_et(k,:) = (-2*pi*k/xl)*et(-k,:)
enddo
end function et_Dx_et
| Function : | |||
| et_Dy_et : | real(8), dimension(-km:km,0:lm)
| ||
| et : | real(8), dimension(-km:km,0:lm), intent(in)
|
入力スペクトルデータに Y 微分(∂y)を作用する.
スペクトルデータの X 微分とは, 対応する格子点データに Y 微分を 作用させたデータのスペクトル変換のことである.
function et_Dy_et(et)
!
! 入力スペクトルデータに Y 微分(∂y)を作用する.
!
! スペクトルデータの X 微分とは, 対応する格子点データに Y 微分を
! 作用させたデータのスペクトル変換のことである.
!
real(8), dimension(-km:km,0:lm) :: et_Dy_et
!(out) スペクトルデータの Y 微分
real(8), dimension(-km:km,0:lm), intent(in) :: et
!(in) 入力スペクトルデータ
et_Dy_et = at_Dy_at(et)
end function et_Dy_et
| Subroutine : | |||
| i : | integer,intent(in)
| ||
| j : | integer,intent(in)
| ||
| k : | integer,intent(in)
| ||
| l : | integer,intent(in)
| ||
| xmin : | real(8),intent(in)
| ||
| xmax : | real(8),intent(in)
| ||
| ymin : | real(8),intent(in)
| ||
| ymax : | real(8),intent(in)
|
スペクトル変換の格子点数, 波数, 領域の大きさを設定する.
他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで 初期設定をしなければならない.
subroutine et_Initial(i,j,k,l,xmin,xmax,ymin,ymax)
!
! スペクトル変換の格子点数, 波数, 領域の大きさを設定する.
!
! 他の関数や変数を呼ぶ前に, 最初にこのサブルーチンを呼んで
! 初期設定をしなければならない.
!
integer,intent(in) :: i ! 格子点数(X)
integer,intent(in) :: j ! 格子点数(Y)
integer,intent(in) :: k ! 切断波数(X)
integer,intent(in) :: l ! 切断波数(Y)
real(8),intent(in) :: xmin, xmax ! X 座標範囲
real(8),intent(in) :: ymin, ymax ! Y 座標範囲
im = i
jm = j
km = k
lm = l
xl = xmax-xmin
yl = ymax-ymin
call ae_initial(im,km,xmin,xmax)
call at_initial(jm,lm,ymin,ymax)
allocate(yx_X(0:jm,0:im-1),yx_Y(0:jm,0:im-1))
yx_X = spread(x_X,1,jm+1)
yx_Y = spread(y_Y,2,im)
end subroutine et_initial
| Function : | |||
| et_Jacobian_et_et : | real(8), dimension(-km:km,0:lm)
| ||
| et_a : | real(8), dimension(-km:km,0:lm), intent(in)
| ||
| et_b : | real(8), dimension(-km:km,0:lm), intent(in)
|
2 つのスペクトルデータからヤコビアン
J(A,B)=(∂xA)(∂yB)-(∂yA)(∂xB)
を計算する.
2 つのスペクトルデータのヤコビアンとは, 対応する 2 つの
格子点データのヤコビアンのスペクトル変換のことである.
function et_Jacobian_et_et(et_a,et_b)
!
! 2 つのスペクトルデータからヤコビアン
!
! J(A,B)=(∂xA)(∂yB)-(∂yA)(∂xB)
!
! を計算する.
!
! 2 つのスペクトルデータのヤコビアンとは, 対応する 2 つの
! 格子点データのヤコビアンのスペクトル変換のことである.
!
real(8), dimension(-km:km,0:lm) :: et_Jacobian_et_et
!(out) 2 つのスペクトルデータのヤコビアン
real(8), dimension(-km:km,0:lm), intent(in) :: et_a
!(in) 1つ目の入力スペクトルデータ
real(8), dimension(-km:km,0:lm), intent(in) :: et_b
!(in) 2つ目の入力スペクトルデータ
et_Jacobian_et_et = et_yx( yx_et(et_Dx_et(et_a)) * yx_et(et_Dy_et(et_b)) -yx_et(et_Dy_et(et_a)) * yx_et(et_Dx_et(et_b)) )
end function et_Jacobian_et_et
| Function : | |||
| et_LaplaInv_et : | real(8), dimension(-km:km,0:lm)
| ||
| et : | real(8), dimension(-km:km,0:lm),intent(in)
| ||
| values : | real(8), dimension(-km:km,2), intent(in), optional
|
境界で一様な値を与える条件(ディリクレ条件)下で, 入力スペクトルデータに逆ラプラシアン(∂xx+∂yy)**(-1)を作用する.
Chebyshev-tau 法による計算
スペクトルデータの逆ラプラシアンとは, 対応する格子点データに 逆ラプラシアンを作用させたデータのスペクトル変換のことである.
function et_LaplaInv_et(et,values)
!
! 境界で一様な値を与える条件(ディリクレ条件)下で,
! 入力スペクトルデータに逆ラプラシアン(∂xx+∂yy)**(-1)を作用する.
!
! Chebyshev-tau 法による計算
!
! スペクトルデータの逆ラプラシアンとは, 対応する格子点データに
! 逆ラプラシアンを作用させたデータのスペクトル変換のことである.
!
real(8), dimension(-km:km,0:lm),intent(in) :: et
!(in) スペクトルデータ
real(8), dimension(-km:km,0:lm) :: et_LaplaInv_et
!(out) スペクトルデータの逆ラプラシアン
real(8), dimension(-km:km,2), intent(in), optional :: values
!(in) 境界値. 省略時は 0 が設定される.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(-km:km,0:lm) :: et_work
real(8), dimension(0:lm,0:lm) :: tt_work
real(8), dimension(0:lm,0:jm) :: ty_work
real(8), dimension(-km:km) :: value1, value2 ! 境界値
logical :: first = .true.
integer :: k,l
save :: alu, kp, first
if (.not. present(values)) then
value1=0
value2=0
else
value1 = values(:,1)
value2 = values(:,2)
endif
if ( first ) then
first = .false.
allocate(alu(-km:km,0:lm,0:lm),kp(-km:km,0:lm))
tt_work=0
do l=0,lm
tt_work(l,l)=1
enddo
ty_work=ay_at(tt_work)
do k=-km,km
alu(k,:,:) = transpose(at_Dy_at(at_Dy_at(tt_work)) - (2*pi*k/xl)**2*tt_work)
alu(k,lm-1,:) = ty_work(:,0)
alu(k,lm,:) = ty_work(:,jm)
enddo
call ludecomp(alu,kp)
endif
et_work = et
et_work(:,lm-1) = value1
et_work(:,lm) = value2
et_LaplaInv_et = lusolve(alu,kp,et_work)
end function et_LaplaInv_et
| Function : | |||
| et_Lapla_et : | real(8), dimension(-km:km,0:lm)
| ||
| et : | real(8), dimension(-km:km,0:lm), intent(in)
|
入力スペクトルデータにラプラシアン(∂xx+∂yy)を作用する.
スペクトルデータのラプラシアンとは, 対応する格子点データに ラプラシアンを作用させたデータのスペクトル変換のことである.
function et_Lapla_et(et)
!
! 入力スペクトルデータにラプラシアン(∂xx+∂yy)を作用する.
!
! スペクトルデータのラプラシアンとは, 対応する格子点データに
! ラプラシアンを作用させたデータのスペクトル変換のことである.
!
real(8), dimension(-km:km,0:lm) :: et_Lapla_et
!(out) スペクトルデータのラプラシアン
real(8), dimension(-km:km,0:lm), intent(in) :: et
!(in) 入力スペクトルデータ
integer k
do k=-km,km
et_Lapla_et(k,:) = -(2*pi*k/xl)**2*et(k,:)
enddo
et_Lapla_et = et_Lapla_et + et_Dy_et(et_Dy_et(et))
end function et_Lapla_et
| Function : | |||
| et_Vor2Strm1_et : | real(8), dimension(-km:km,0:lm)
| ||
| et : | real(8), dimension(-km:km,0:lm),intent(in)
| ||
| values : | real(8), dimension(2), intent(in), optional
| ||
| rigid : | logical, dimension(2), intent(in), optional |
渦度から流線を求める.
注意 : 以下の点がうまくなく現在保留中. 使用禁止
* Y 方向の切断波数を lm=jm と設定しておかないと
計算結果が振動して不安定となる.
* 非圧縮流体計算の拡散問題の時間発展が安定に計算できない.
Chebyshev-tauとこの問題の相性は良くないらしい.
Chebyshev-tau 法による計算 渦度 \zeta を与えて流線 \psi を求める.
\nabla^2 \psi = \zeta, \psi = const. at boundaries.
粘着条件
\DP{\psi}{y} = 0 at boundaries
応力なし条件
\DP[2]{\psi}{y} = 0 at boundaries
デフォルトは粘着条件
\nabla^4 \psi = \nabla^2\zeta を解く 4 階の微分方程式にすることで境界条件の数とつじつまがあう.
function et_Vor2Strm1_et(et,values,rigid)
! 渦度から流線を求める.
!
! 注意 : 以下の点がうまくなく現在保留中. 使用禁止
! * Y 方向の切断波数を lm=jm と設定しておかないと
! 計算結果が振動して不安定となる.
! * 非圧縮流体計算の拡散問題の時間発展が安定に計算できない.
! Chebyshev-tauとこの問題の相性は良くないらしい.
!
! Chebyshev-tau 法による計算
! 渦度 \zeta を与えて流線 \psi を求める.
! \nabla^2 \psi = \zeta,
! \psi = const. at boundaries.
! 粘着条件
! \DP{\psi}{y} = 0 at boundaries
! 応力なし条件
! \DP[2]{\psi}{y} = 0 at boundaries
!
! デフォルトは粘着条件
!
! \nabla^4 \psi = \nabla^2\zeta を解く
! 4 階の微分方程式にすることで境界条件の数とつじつまがあう.
!
!
real(8), dimension(-km:km,0:lm),intent(in) :: et
!(in) 渦度
real(8), dimension(-km:km,0:lm) :: et_Vor2Strm1_et
!(out) 流線境界値. 境界で一定なので波数 0 成分のみ
real(8), dimension(2), intent(in), optional :: values
!(in) 境界条件スイッチ
logical, dimension(2), intent(in), optional :: rigid
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(-km:km,0:lm) :: et_work
real(8), dimension(0:lm,0:lm) :: tt_work
real(8), dimension(0:lm,0:jm) :: ty_work
real(8) :: value1, value2 ! 境界値
logical :: rigid1, rigid2 ! 境界条件
logical :: first = .true.
integer :: k,l
save :: alu, kp, first
if (.not. present(values)) then
value1=0
value2=0
else
value1 = values(1)
value2 = values(2)
endif
if (.not. present(rigid)) then
rigid1=.true.
rigid2=.true.
else
rigid1 = rigid(1)
rigid2 = rigid(2)
endif
if ( first ) then
first = .false.
allocate(alu(-km:km,0:lm,0:lm),kp(-km:km,0:lm))
tt_work=0
do l=0,lm
tt_work(l,l)=1
enddo
do k=-km,km
alu(k,:,:) = transpose( at_Dy_at(at_Dy_at(at_Dy_at(at_Dy_at(tt_work)))) - 2 * (2*pi*k/xl)**2 * at_Dy_at(at_Dy_at(tt_work)) + (2*pi*k/xl)**4*tt_work )
enddo
! 運動学的条件. 流線は境界で一定
ty_work=ay_at(tt_work)
do k=-km,km
alu(k,lm-1,:) = ty_work(:,0)
alu(k,lm,:) = ty_work(:,jm)
enddo
! 力学的条件粘着条件
if ( rigid1 ) then
ty_work=ay_at(at_Dy_at(tt_work))
else
ty_work=ay_at(at_Dy_at(at_Dy_at(tt_work)))
endif
do k=-km,km
alu(k,lm-3,:) = ty_work(:,0)
enddo
! 力学的条件粘着条件
if ( rigid2 ) then
ty_work=ay_at(at_Dy_at(tt_work))
else
ty_work=ay_at(at_Dy_at(at_Dy_at(tt_work)))
endif
do k=-km,km
alu(k,lm-2,:) = ty_work(:,jm)
enddo
call ludecomp(alu,kp)
endif
et_work = et_Lapla_et(et)
et_work(:,lm-3) = 0 ! 力学的条件
et_work(:,lm-2) = 0 ! 力学的条件
et_work(:,lm-1) = 0 ! 運動学的条件. 波数 0 以外は 0
et_work(0,lm-1) = value1*2 ! 運動学的条件. 波数 0 は重み 1/2
et_work(:,lm) = 0 ! 運動学的条件. 波数 0 以外は 0
et_work(0,lm) = value2*2 ! 運動学的条件. 波数 0 は重み 1/2
et_Vor2Strm1_et = lusolve(alu,kp,et_work)
end function et_Vor2Strm1_et
| Function : | |||
| et_Vor2Strm_et : | real(8), dimension(-km:km,0:lm) | ||
| et : | real(8), dimension(-km:km,0:lm),intent(in) | ||
| values : | real(8), dimension(2), intent(in), optional
| ||
| rigid : | logical, dimension(2), intent(in), optional |
渦度から流線を求める.
注意 : 以下の点がうまくなく現在保留中. 使用禁止
* Y 方向の切断波数を lm=jm と設定しておかないと
計算結果が振動して不安定となる.
* 非圧縮流体計算の拡散問題の時間発展が安定に計算できない.
Chebyshev-tauとこの問題の相性は良くないらしい.
Chebyshev-tau 法による計算 渦度 \zeta を与えて流線 \psi を求める.
\nabla^2 \psi = \zeta, \psi = const. at boundaries.
粘着条件
\DP{\psi}{y} = 0 at boundaries
応力なし条件
\DP[2]{\psi}{y} = 0 at boundaries
l=0,1,lm-1,lm 成分の式の代わりに境界条件を与える. 渦度の低次成分を無視することは \nabla^4 \psi = \zeta^2 を解いていることに相当. 4 階の微分方程式にすることで境界条件の数とつじつまがあう.
function et_Vor2Strm_et(et,values,rigid)
!
! 渦度から流線を求める.
!
! 注意 : 以下の点がうまくなく現在保留中. 使用禁止
! * Y 方向の切断波数を lm=jm と設定しておかないと
! 計算結果が振動して不安定となる.
! * 非圧縮流体計算の拡散問題の時間発展が安定に計算できない.
! Chebyshev-tauとこの問題の相性は良くないらしい.
!
! Chebyshev-tau 法による計算
! 渦度 \zeta を与えて流線 \psi を求める.
! \nabla^2 \psi = \zeta,
! \psi = const. at boundaries.
! 粘着条件
! \DP{\psi}{y} = 0 at boundaries
! 応力なし条件
! \DP[2]{\psi}{y} = 0 at boundaries
!
! l=0,1,lm-1,lm 成分の式の代わりに境界条件を与える.
! 渦度の低次成分を無視することは
! \nabla^4 \psi = \zeta^2 を解いていることに相当.
! 4 階の微分方程式にすることで境界条件の数とつじつまがあう.
real(8), dimension(-km:km,0:lm),intent(in) :: et
real(8), dimension(-km:km,0:lm) :: et_Vor2Strm_et
! 流線境界値. 境界で一定なので波数 0 成分のみ
real(8), dimension(2), intent(in), optional :: values
! 境界条件スイッチ
logical, dimension(2), intent(in), optional :: rigid
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(-km:km,0:lm) :: et_work
real(8), dimension(0:lm,0:lm) :: tt_work
real(8), dimension(0:lm,0:jm) :: ty_work
real(8) :: value1, value2 ! 境界値
logical :: rigid1, rigid2 ! 境界条件
logical :: first = .true.
integer :: k,l
save :: alu, kp, first
if (.not. present(values)) then
value1=0
value2=0
else
value1 = values(1)
value2 = values(2)
endif
if (.not. present(rigid)) then
rigid1=.true.
rigid2=.true.
else
rigid1 = rigid(1)
rigid2 = rigid(2)
endif
if ( first ) then
first = .false.
allocate(alu(-km:km,0:lm,0:lm),kp(-km:km,0:lm))
tt_work=0
do l=0,lm
tt_work(l,l)=1
enddo
do k=-km,km
alu(k,:,:) = transpose(at_Dy_at(at_Dy_at(tt_work)) - (2*pi*k/xl)**2*tt_work)
enddo
! 運動学的条件. 流線は境界で一定
ty_work=ay_at(tt_work)
do k=-km,km
alu(k,lm-1,:) = ty_work(:,0)
alu(k,lm,:) = ty_work(:,jm)
enddo
! 力学的条件粘着条件
if ( rigid1 ) then
ty_work=ay_at(at_Dy_at(tt_work))
else
ty_work=ay_at(at_Dy_at(at_Dy_at(tt_work)))
endif
do k=-km,km
alu(k,0,:) = ty_work(:,0)
enddo
! 力学的条件粘着条件
if ( rigid2 ) then
ty_work=ay_at(at_Dy_at(tt_work))
else
ty_work=ay_at(at_Dy_at(at_Dy_at(tt_work)))
endif
do k=-km,km
alu(k,1,:) = ty_work(:,jm)
enddo
call ludecomp(alu,kp)
endif
et_work = et
et_work(:,0) = 0 ! 力学的条件
et_work(:,1) = 0 ! 力学的条件
et_work(:,lm-1) = 0 ! 運動学的条件. 波数 0 以外は 0
et_work(0,lm-1) = value1*2 ! 運動学的条件. 波数 0 は重み 1/2
et_work(:,lm) = 0 ! 運動学的条件. 波数 0 以外は 0
et_work(0,lm) = value2*2 ! 運動学的条件. 波数 0 は重み 1/2
et_Vor2Strm_et = lusolve(alu,kp,et_work)
end function et_Vor2Strm_et
| Function : | |||
| et_yx : | real(8), dimension(-km:km,0:lm)
| ||
| yx : | real(8), dimension(0:jm,0:im-1), intent(in)
|
格子データからスペクトルデータへ変換する.
function et_yx(yx)
!
! 格子データからスペクトルデータへ変換する.
!
real(8), dimension(-km:km,0:lm) :: et_yx
!(out) スペクトルデータ
real(8), dimension(0:jm,0:im-1), intent(in) :: yx
!(in) 格子点データ
et_yx = at_ay(transpose(ae_ax(yx)))
end function et_yx
| Function : | |||
| ey_Vor2Strm_ey : | real(8), dimension(-km:km,0:jm)
| ||
| ey : | real(8), dimension(-km:km,0:jm),intent(in)
| ||
| values : | real(8), dimension(2), intent(in), optional
| ||
| cond : | character(len=2), intent(in), optional
|
渦度から流線を求める.
渦度から流線を求める. Y 方向格子点空間での Chebyshev-Collocation 法による計算
渦度 \zeta を与えて流線 \psi を求める.
\nabla^2 \psi = \zeta, \psi = const. at boundaries.
粘着条件
\DP{\psi}{y} = 0 at boundaries
応力なし条件
\DP[2]{\psi}{y} = 0 at boundaries
注 : 最初に呼ばれるときの境界条件で以後計算される(要仕様変更)
function ey_Vor2Strm_ey(ey,values,cond) ! 渦度から流線を求める.
!
! 渦度から流線を求める.
! Y 方向格子点空間での Chebyshev-Collocation 法による計算
!
! 渦度 \zeta を与えて流線 \psi を求める.
! \nabla^2 \psi = \zeta,
! \psi = const. at boundaries.
! 粘着条件
! \DP{\psi}{y} = 0 at boundaries
! 応力なし条件
! \DP[2]{\psi}{y} = 0 at boundaries
!
! 注 : 最初に呼ばれるときの境界条件で以後計算される(要仕様変更)
real(8), dimension(-km:km,0:jm),intent(in) :: ey
!(in) 入力渦度分布
real(8), dimension(-km:km,0:jm) :: ey_Vor2Strm_ey
!(out) 出力流線分布
real(8), dimension(2), intent(in), optional :: values
!(in) 流線境界値. 境界で一定なので波数 0 成分のみ
! 省略時は 0.
character(len=2), intent(in), optional :: cond
! (in) 境界条件スイッチ. 省略時は 'RR'
! RR : 両端粘着条件
! RF,FR : 粘着/応力なし条件
! FF : 両端応力なし条件
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(-km:km,0:jm) :: ey_work
real(8), dimension(0:jm,0:jm) :: yy
real(8), dimension(0:jm,0:jm) :: yy_work
real(8) :: value1, value2 ! 境界値
logical :: rigid1, rigid2 ! 境界条件
logical :: first = .true.
integer :: k,j
save :: alu, kp, first
if (.not. present(values)) then
value1=0
value2=0
else
value1 = values(1)
value2 = values(2)
endif
if (.not. present(cond)) then
rigid1=.TRUE.
rigid2=.TRUE.
else
select case (cond)
case ('RR')
rigid1 = .TRUE.
rigid2 = .TRUE.
case ('RF')
rigid1 = .TRUE.
rigid2 = .FALSE.
case ('FR')
rigid1 = .FALSE.
rigid2 = .TRUE.
case ('FF')
rigid1 = .FALSE.
rigid2 = .FALSE.
case default
call MessageNotify('E','ey_Vor2Strm_ey','B.C. not supported')
end select
endif
if ( first ) then
first = .false.
allocate(alu(-km:km,0:jm,0:jm),kp(-km:km,0:jm))
yy=0
do j=0,jm
yy(j,j)=1
enddo
do k=-km,km
alu(k,:,:) = transpose( ay_at(at_Dy_at(at_Dy_at(at_ay(yy)))) - (2*pi*k/xl)** 2* yy )
enddo
! 運動学的条件. 流線は境界で一定
yy_work=yy
do k=-km,km
alu(k,0,:) = yy_work(:,0)
alu(k,jm,:) = yy_work(:,jm)
enddo
! 力学的条件粘着条件
if ( rigid1 ) then
yy_work=ay_at(at_Dy_at(at_ay(yy)))
else
yy_work=ay_at(at_Dy_at(at_Dy_at(at_ay(yy))))
endif
do k=-km,km
alu(k,1,:) = yy_work(:,0)
enddo
! 力学的条件粘着条件
if ( rigid2 ) then
yy_work=ay_at(at_Dy_at(at_ay(yy)))
else
yy_work=ay_at(at_Dy_at(at_Dy_at(at_ay(yy))))
endif
do k=-km,km
alu(k,jm-1,:) = yy_work(:,jm)
enddo
call ludecomp(alu,kp)
endif
ey_work = ey
ey_work(:,1) = 0 ! 力学的条件
ey_work(:,jm-1) = 0 ! 力学的条件
ey_work(:,0) = 0 ! 運動学的条件. 波数 0 以外は 0
ey_work(0,0) = value1*2 ! 運動学的条件. 波数 0 は重み 1/2
ey_work(:,jm) = 0 ! 運動学的条件. 波数 0 以外は 0
ey_work(0,jm) = value2*2 ! 運動学的条件. 波数 0 は重み 1/2
ey_Vor2Strm_ey = lusolve(alu,kp,ey_work)
end function ey_Vor2Strm_ey
| Function : | |||
| t_Dx_t : | real(8), dimension(size(t_data))
| ||
| t_data : | real(8), dimension(:), intent(in)
|
入力チェビシェフデータに X 微分を作用する(1 次元配列用).
チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.
The entity is at_module#t_Dx_t
| Function : | |||
| t_g : | double precision, dimension(0:km)
| ||
| g_data : | double precision, dimension(:), intent(in)
|
台形格子 -> スペクトル
格子データからチェビシェフデータへ変換する(2 次元配列用).
The entity is at_module#t_g
| Function : | |||
| x_AvrY_yx : | real(8), dimension(0:im-1)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
2 次元格子点データの Y 方向平均
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.
function x_AvrY_yx(yx)
!
! 2 次元格子点データの Y 方向平均
!
! 実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し,
! y_Y_Weight の総和で割ることで平均している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8), dimension(0:im-1) :: x_AvrY_yx
!(out) 平均された 1 次元(X)格子点
x_AvrY_yx = x_IntY_yx(yx)/sum(y_Y_weight)
end function x_AvrY_yx
| Function : | |||
| x_IntY_yx : | real(8), dimension(0:im-1)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
Y 方向積分
2 次元格子点データの Y 方向積分
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
function x_IntY_yx(yx) ! Y 方向積分
!
! 2 次元格子点データの Y 方向積分
!
! 実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8), dimension(0:im-1) :: x_IntY_yx
!(out) 積分された 1 次元(X)格子点データ
integer :: j
! 作業変数
x_IntY_yx = 0.0d0
do j=0,jm
x_IntY_yx(:) = x_IntY_yx(:) + yx(j,:) * y_Y_Weight(j)
enddo
end function x_IntY_yx
| Variable : | |||
| g_X(:) : | real(8), allocatable
|
The entity is at_module#g_X
| Variable : | |||
| g_X_Weight(:) : | real(8), allocatable
|
The entity is at_module#g_X_Weight
| Variable : | |||
| g_x_weight(:) : | real(8), allocatable
|
The entity is ae_module#g_x_weight
| Function : | |||
| g_e : | real(8), dimension(0:im-1)
| ||
| e : | real(8), dimension(-km:km), intent(in)
|
スペクトルデータから格子点データへ変換する(1 次元データ用)
The entity is ae_module#g_e
| Function : | |||
| y_AvrX_yx : | real(8), dimension(0:jm)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
2 次元格子点データの X 方向平均
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.
function y_AvrX_yx(yx)
!
! 2 次元格子点データの X 方向平均
!
! 実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し,
! x_X_Weight の総和で割ることで平均している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8), dimension(0:jm) :: y_AvrX_yx
!(out) 平均された 1 次元(Y)格子点
y_AvrX_yx = y_IntX_yx(yx)/sum(x_X_weight)
end function y_AvrX_yx
| Function : | |||
| y_IntX_yx : | real(8), dimension(0:jm)
| ||
| yx : | real(8), dimension(0:jm,0:im-1)
|
X 方向積分
2 次元格子点データの X 方向積分
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
function y_IntX_yx(yx) ! X 方向積分
!
! 2 次元格子点データの X 方向積分
!
! 実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
!
real(8), dimension(0:jm,0:im-1) :: yx
!(in) 2 次元格子点データ
real(8), dimension(0:jm) :: y_IntX_yx
!(out) 積分された 1 次元(Y)格子点データ
integer :: i
! 作業変数
y_IntX_yx = 0.0d0
do i=0,im-1
y_IntX_yx(:) = y_IntX_yx(:) + yx(:,i) * x_X_Weight(i)
enddo
end function y_IntX_yx
| Variable : | |||
| g_X(:) : | real(8), allocatable
|
The entity is at_module#g_X
| Variable : | |||
| g_X_Weight(:) : | real(8), allocatable
|
The entity is at_module#g_X_Weight
| Variable : | |||
| g_x_weight(:) : | real(8), allocatable
|
The entity is ae_module#g_x_weight
| Function : | |||
| g_t : | double precision, dimension(0:im)
| ||
| t_data : | double precision, dimension(:), intent(in)
|
チェビシェフデータから格子データへ変換する(1 次元配列用).
The entity is at_module#g_t
| Function : | |||
| yx_et : | real(8), dimension(0:jm,0:im-1)
| ||
| et : | real(8), dimension(-km:km,0:lm), intent(in)
|
スペクトルデータから格子データへ変換する.
function yx_et(et)
!
! スペクトルデータから格子データへ変換する.
!
real(8), dimension(0:jm,0:im-1) :: yx_et
!(out) 格子点データ
real(8), dimension(-km:km,0:lm), intent(in) :: et
!(in) スペクトルデータ
yx_et = ax_ae(transpose(ay_at(et)))
end function yx_et