Class | xyz_deriv_c4_module |
In: |
../src/utils/xyz_deriv_c4_module.f90
|
xyz_deriv_module は, 3 次元 (xyz 方向) 不等間隔交互格子を用いた有限差分 法に基づく数値モデルのための, 微分演算を行う Fortran90 関数を提供する. 微分演算は 4 次精度中心差分を用いて行う.
このモジュールは下請けモジュールとして xyz_base_module, data_type モジュールを用いている.
モジュール内の変数と手続きの命名方法については xyz_module を参照のこと.
xyz_dx_pyz, pyz_dx_xyz : | x 方向 1 階微分を計算する |
xyz_dy_xqz, xqz_dy_xyz : | y 方向 1 階微分を計算する |
xyz_dz_xyr, xyr_dz_xyz : | z 方向 1 階微分を計算する |
差分式は格子間隔を一定とした
du/dx(i) = 9/8*[u(i)-u(i-1)]/dx(i) - 1/24*[u(i+1)-u(i-2)]/dx(i)
を用いている. 不等間隔格子の場合は精度が落ちることに注意.
Function : | |
paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
xaa_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function paa_c4dx_xaa(xaa_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: xaa_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) integer :: ix ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! paa_c4dx_xaa = 0.0d0 ! 1 階微分の計算 ! do ix = imin+1, imax-2 paa_c4dx_xaa(ix,:,:) = ((xaa_Var(ix+1,:,:) - xaa_Var(ix,:,:))*9.0d0/(8.0d0*p_dx(ix)) - ((xaa_Var(ix+2,:,:) - xaa_Var(ix-1,:,:))/(24.0d0*p_dx(ix)))) end do end function paa_c4dx_xaa
Function : | |
aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aya_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aqa_c4dy_aya(aya_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aya_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) integer :: jy ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aqa_c4dy_aya = 0.0d0 ! y 方向一様な場合は微分はゼロ. if (jmin == jmax) return ! 1 階微分の計算 ! do jy = jmin+1, jmax-2 aqa_c4dy_aya(:,jy,:) = ((aya_Var(:,jy+1,:) - aya_Var(:,jy,:))*9.0d0/(8.0d0*q_dy(jy)) - ((aya_Var(:,jy+2,:) - aya_Var(:,jy-1,:))/(24.0d0*q_dy(jy)))) end do end function aqa_c4dy_aya
Function : | |
paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
xaa_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function paa_c4dx_xaa(xaa_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: xaa_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) integer :: ix ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! paa_c4dx_xaa = 0.0d0 ! 1 階微分の計算 ! do ix = imin+1, imax-2 paa_c4dx_xaa(ix,:,:) = ((xaa_Var(ix+1,:,:) - xaa_Var(ix,:,:))*9.0d0/(8.0d0*p_dx(ix)) - ((xaa_Var(ix+2,:,:) - xaa_Var(ix-1,:,:))/(24.0d0*p_dx(ix)))) end do end function paa_c4dx_xaa
Function : | |
aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aaz_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aar_c4dz_aaz(aaz_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aaz_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) integer :: kz ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aar_c4dz_aaz = 0.0d0 ! 1 階微分の計算 ! * r_c4dz_z を用いて計算 ! do kz = kmin+1, kmax-2 aar_c4dz_aaz(:,:,kz) = ((aaz_Var(:,:,kz+1) - aaz_Var(:,:,kz))*9.0d0/(8.0d0*r_dz(kz)) - ((aaz_Var(:,:,kz+2) - aaz_Var(:,:,kz-1))/(24.0d0*r_dz(kz)))) end do end function aar_c4dz_aaz
Function : | |
paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
xaa_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function paa_c4dx_xaa(xaa_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: xaa_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: paa_c4dx_xaa(imin:imax,jmin:jmax,kmin:kmax) integer :: ix ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! paa_c4dx_xaa = 0.0d0 ! 1 階微分の計算 ! do ix = imin+1, imax-2 paa_c4dx_xaa(ix,:,:) = ((xaa_Var(ix+1,:,:) - xaa_Var(ix,:,:))*9.0d0/(8.0d0*p_dx(ix)) - ((xaa_Var(ix+2,:,:) - xaa_Var(ix-1,:,:))/(24.0d0*p_dx(ix)))) end do end function paa_c4dx_xaa
Function : | |
aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aya_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aqa_c4dy_aya(aya_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aya_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) integer :: jy ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aqa_c4dy_aya = 0.0d0 ! y 方向一様な場合は微分はゼロ. if (jmin == jmax) return ! 1 階微分の計算 ! do jy = jmin+1, jmax-2 aqa_c4dy_aya(:,jy,:) = ((aya_Var(:,jy+1,:) - aya_Var(:,jy,:))*9.0d0/(8.0d0*q_dy(jy)) - ((aya_Var(:,jy+2,:) - aya_Var(:,jy-1,:))/(24.0d0*q_dy(jy)))) end do end function aqa_c4dy_aya
Function : | |
aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aaz_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aar_c4dz_aaz(aaz_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aaz_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) integer :: kz ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aar_c4dz_aaz = 0.0d0 ! 1 階微分の計算 ! * r_c4dz_z を用いて計算 ! do kz = kmin+1, kmax-2 aar_c4dz_aaz(:,:,kz) = ((aaz_Var(:,:,kz+1) - aaz_Var(:,:,kz))*9.0d0/(8.0d0*r_dz(kz)) - ((aaz_Var(:,:,kz+2) - aaz_Var(:,:,kz-1))/(24.0d0*r_dz(kz)))) end do end function aar_c4dz_aaz
Function : | |
aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aya_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aqa_c4dy_aya(aya_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aya_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aqa_c4dy_aya(imin:imax,jmin:jmax,kmin:kmax) integer :: jy ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aqa_c4dy_aya = 0.0d0 ! y 方向一様な場合は微分はゼロ. if (jmin == jmax) return ! 1 階微分の計算 ! do jy = jmin+1, jmax-2 aqa_c4dy_aya(:,jy,:) = ((aya_Var(:,jy+1,:) - aya_Var(:,jy,:))*9.0d0/(8.0d0*q_dy(jy)) - ((aya_Var(:,jy+2,:) - aya_Var(:,jy-1,:))/(24.0d0*q_dy(jy)))) end do end function aqa_c4dy_aya
Function : | |
aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
aaz_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function aar_c4dz_aaz(aaz_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: aaz_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: aar_c4dz_aaz(imin:imax,jmin:jmax,kmin:kmax) integer :: kz ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! aar_c4dz_aaz = 0.0d0 ! 1 階微分の計算 ! * r_c4dz_z を用いて計算 ! do kz = kmin+1, kmax-2 aar_c4dz_aaz(:,:,kz) = ((aaz_Var(:,:,kz+1) - aaz_Var(:,:,kz))*9.0d0/(8.0d0*r_dz(kz)) - ((aaz_Var(:,:,kz+2) - aaz_Var(:,:,kz-1))/(24.0d0*r_dz(kz)))) end do end function aar_c4dz_aaz
Function : | |
xyz_c4dx_pyz(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
pyz_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function xyz_c4dx_pyz(pyz_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: pyz_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: xyz_c4dx_pyz(imin:imax,jmin:jmax,kmin:kmax) integer :: ix ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! xyz_c4dx_pyz = 0.0d0 ! 1 階微分の計算 ! do ix = imin+2, imax-1 xyz_c4dx_pyz(ix,:,:) = ((pyz_Var(ix,:,:) - pyz_Var(ix-1,:,:))*9.0d0/(8.0d0*x_dx(ix)) - ((pyz_Var(ix+1,:,:) - pyz_Var(ix-2,:,:))/(24.0d0*x_dx(ix)))) end do end function xyz_c4dx_pyz
Function : | |
xyz_c4dy_xqz(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
xqz_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function xyz_c4dy_xqz(xqz_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: xqz_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: xyz_c4dy_xqz(imin:imax,jmin:jmax,kmin:kmax) integer :: jy ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! xyz_c4dy_xqz = 0.0d0 ! y 方向一様な場合は微分はゼロ. if (jmin == jmax) return ! 1 階微分の計算 ! do jy = jmin+2, jmax-1 xyz_c4dy_xqz(:,jy,:) = ((xqz_Var(:,jy,:) - xqz_Var(:,jy-1,:))*9.0d0/(8.0d0*y_dy(jy)) - ((xqz_Var(:,jy+1,:) - xqz_Var(:,jy-2,:))/(24.0d0*y_dy(jy)))) end do end function xyz_c4dy_xqz
Function : | |
xyz_c4dz_xyr(imin:imax,jmin:jmax,kmin:kmax) : | real(DP) |
xyr_Var(imin:imax,jmin:jmax,kmin:kmax) : | real(DP),intent(in) |
半整数格子点上の 1 階微分を計算する
function xyz_c4dz_xyr(xyr_Var) ! 半整数格子点上の 1 階微分を計算する real(DP),intent(in) :: xyr_Var(imin:imax,jmin:jmax,kmin:kmax) real(DP) :: xyz_c4dz_xyr(imin:imax,jmin:jmax,kmin:kmax) integer :: kz ! 初期化 ! * 0 割りを防ぐためにはマシンイプシロン値を用いるべき ! xyz_c4dz_xyr = 0.0d0 ! 1 階微分の計算 ! do kz = kmin+2, kmax-1 xyz_c4dz_xyr(:,:,kz) = ((xyr_Var(:,:,kz) - xyr_Var(:,:,kz-1))*9.0d0/(8.0d0*z_dz(kz)) - ((xyr_Var(:,:,kz+1) - xyr_Var(:,:,kz-2))/(24.0d0*z_dz(kz)))) end do end function xyz_c4dz_xyr