Class | typhoon_analy |
In: |
typhoon_analy.f90
|
台風感度実験用スペシャル解析モジュール
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
fg(2) : | integer, intent(in)
| ||
pres(size(x),size(y)) : | real, intent(in)
| ||
search_dis : | real, intent(in)
| ||
var_dis : | real, intent(in)
| ||
center(2) : | integer, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
stdopt : | logical, intent(in), optional
|
Braun (2002) の方法を基に台風の中心を推定する.
subroutine DC_Braun( x, y, fg, pres, search_dis, var_dis, center, undef, stdopt ) ! Braun (2002) の方法を基に台風の中心を推定する. use Math_Const use Max_Min use Statistics implicit none real, intent(in) :: x(:) ! x 方向の座標 [m] real, intent(in) :: y(:) ! y 方向の座標 [m] integer, intent(in) :: fg(2) ! 中心点の第一推定値 (通常気圧の最低値等から得る) ! fg(1) = x 方向の要素番号, fg(2) = y 方向の要素番号 real, intent(in) :: pres(size(x),size(y)) ! ある高度での気圧 (地表面気圧でもよい.) ! ただし, 地表面気圧の場合は, 海面校正しておくこと. real, intent(in) :: search_dis ! 検索する領域 (fg の位置を中心に) ! 例えば, 100000.0 なら, fg を中心に縦横 100 km 四方 real, intent(in) :: var_dis ! 推定中心位置から偏差を計算する半径 [m] integer, intent(inout) :: center(2) ! 求めた中心点の各要素数 real, intent(in), optional :: undef ! 気圧に未定義値がある場合, その未定義値. ! 本ルーチンでは, 未定義値がある場合, ! その格子点のみ偏差計算に使用しない. logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) integer :: i, j, ix, jy, nx, ny, nr, nt, nxgmin, nxgmax, nygmin, nygmax integer :: ompnum, tmp_o_num real :: dr, dtheta, undeff, tmpmin, tmp_anom, tmp_counter real, allocatable, dimension(:) :: rad, theta real, allocatable, dimension(:,:) :: anom_check real, allocatable, dimension(:,:,:) :: apres logical :: stderr !-- OpenMP 用整数関数 !$ integer :: OMP_GET_THREAD_NUM, OMP_GET_MAX_THREADS nx=size(x) ny=size(y) if(present(undef))then undeff=undef else undeff=-999.0 end if if(present(stdopt))then stderr=stdopt else stderr=.false. end if !-- 検索格子点範囲を規程 call interpo_search_2d( x, y, x(fg(1))-0.5*search_dis, y(fg(2))-0.5*search_dis, nxgmin, nygmin, undeff=int(undeff), stdopt=stderr ) call interpo_search_2d( x, y, x(fg(1))+0.5*search_dis, y(fg(2))+0.5*search_dis, nxgmax, nygmax, undeff=int(undeff), stdopt=stderr ) if(nxgmin==int(undeff))then ! 領域外参照の場合の処理 nxgmin=1 end if if(nygmin==int(undeff))then ! 領域外参照の場合の処理 nygmin=1 end if if(nxgmax==int(undeff))then ! 領域外参照の場合の処理 nxgmax=1 end if if(nygmax==int(undeff))then ! 領域外参照の場合の処理 nygmax=1 end if ! nxgmin=fg(1)-(search_dis-1)/2 ! nxgmax=fg(1)+(search_dis-1)/2 ! nygmin=fg(2)-(search_dis-1)/2 ! nygmax=fg(2)+(search_dis-1)/2 allocate(anom_check(nx,ny)) !-- openmp での条件付きコンパイル !-- 接線平均アノマリの箇所を openmp 並列したいが, !-- apres が inout 属性なので, private 属性を指定しないと !-- thread ごとに apres の値が上書きされてしまう. !-- そこで, threads number を参照した 3 次元配列にして, !-- thread ごとに別の配列を使うように変更. ompnum=1 !$ ompnum=OMP_GET_MAX_THREADS() ! OpenMP が有効の場合はここも有効. allocate(apres(nx,ny,ompnum)) !-- 円筒系への変換の際には, var_dis での接線解像度が x, y に等しくなるように !-- 設定する. nr=int(var_dis/(x(2)-x(1)))+1 nt=int(2.0*pi*var_dis/(x(2)-x(1)))+1 ! nt=4 dr=x(2)-x(1) dtheta=2.0*pi/real(nt-1) allocate(rad(nr)) allocate(theta(nt)) rad=(/((dr*real(i-1)),i=1,nr)/) theta=(/((dtheta*real(i-1)),i=1,nt)/) !-- 各探索格子点について, 接線平均偏差をとり, !-- 各格子点における偏差の合計を計算する. anom_check=undeff ! 探索範囲外の格子点にはすべて undeff を入れる. ! 探索範囲内の格子にはあとでゼロが入れられ初期化する. tmp_o_num=1 ! OpenMP が有効でない場合, この値が apres の 3 次元目へ !$omp parallel default(shared) !$omp do schedule(dynamic) private(i,j,tmp_anom,tmp_o_num,tmp_counter) do j=nygmin,nygmax do i=nxgmin,nxgmax if(pres(i,j)/=undeff)then !$ tmp_o_num=OMP_GET_THREAD_NUM()+1 ! OpenMP が有効の場合, この値が apres の 3 次元目へ apres(:,:,tmp_o_num)=0.0 tmp_counter=0 call tangent_mean_anom_scal_Cart( x, y, x(i), y(j), pres, rad, theta, apres(:,:,tmp_o_num), undef=undeff, undefg=undeff, undefgc='inc', stdopt=stderr ) tmp_anom=0.0 do jy=1,ny do ix=1,nx if(apres(ix,jy,tmp_o_num)/=undeff)then tmp_counter=tmp_counter+1 tmp_anom=tmp_anom+apres(ix,jy,tmp_o_num)*apres(ix,jy,tmp_o_num) end if end do end do if(tmp_counter>0)then ! 平均したが undef しかないときは更新しない. anom_check(i,j)=tmp_anom end if end if end do end do !$omp end do !$omp end parallel !-- 計算した偏差の合計値のうち, 最小となる格子点を求める. call min_val_2d( anom_check, center(1), center(2), tmpmin, undef=undeff ) if(stderr.eqv..false.)then if(center(1)==nx+1.or.center(2)==ny+1)then write(*,*) "*** WARNING *** : DC_Braun (typhoon_analy)" write(*,*) "Setting the undef only point." end if end if deallocate(rad) deallocate(theta) deallocate(anom_check) deallocate(apres) end subroutine DC_Braun
Subroutine : | |||
r(:) : | real, intent(in)
| ||
coril(size(r)) : | real, intent(in)
| ||
v(size(r)) : | real, intent(in)
| ||
rho(size(r)) : | real, intent(in)
| ||
r_ref : | real, intent(in)
| ||
p_ref : | real, intent(in)
| ||
pres(size(r)) : | real, intent(inout)
|
傾度風平衡場を満たす気圧場を計算する.
subroutine grad_wind_pres( r, coril, v, rho, r_ref, p_ref, pres ) ! 傾度風平衡場を満たす気圧場を計算する. use Algebra implicit none real, intent(in) :: r(:) ! r 方向の位置座標 [m] real, intent(in) :: coril(size(r)) ! コリオリパラメータ [/s] real, intent(in) :: v(size(r)) ! r 方向の位置座標 [m] real, intent(in) :: rho(size(r)) ! 密度 [kg/m^3] real, intent(in) :: r_ref ! 積分定数となる位置座標 [m] real, intent(in) :: p_ref ! r_ref での気圧 (積分定数) [Pa] real, intent(inout) :: pres(size(r)) ! 傾度風平衡での気圧 [Pa] integer :: i, nr real :: grad(size(r)) nr=size(r) do i=1,nr if(r(i)/=0.0)then grad(i)=rho(i)*(v(i)*v(i)/r(i)+coril(i)*v(i)) else grad(i)=0.0 end if end do do i=1,nr if(r(i)<r_ref)then call rectangle_int( r, grad, r(i), r_ref, pres(i) ) pres(i)=p_ref-pres(i) else if(r(i)>r_ref)then call rectangle_int( r, grad, r_ref, r(i), pres(i) ) pres(i)=p_ref+pres(i) else pres(i)=p_ref end if end if end do end subroutine
Subroutine : | |||
r(:) : | real, intent(in)
| ||
z(:) : | real, intent(in)
| ||
coril(size(r),size(z)) : | real, intent(in)
| ||
v(size(r),size(z)) : | real, intent(in)
| ||
pres_s(size(z)) : | real, intent(in)
| ||
rho_s(size(z)) : | real, intent(in)
| ||
pres(size(r),size(z)) : | real, intent(inout)
| ||
rho(size(r),size(z)) : | real, intent(inout)
| ||
error : | real, intent(in), optional
| ||
dl : | integer, intent(in), optional
|
サウンディングと軸対称流から静力学・傾度風平衡場の計算.
subroutine hydro_grad_eqb( r, z, coril, v, pres_s, rho_s, pres, rho, error, dl ) ! サウンディングと軸対称流から静力学・傾度風平衡場の計算. use Thermo_Const use Phys_Const use algebra use stdio use Derivation implicit none real, intent(in) :: r(:) ! 動径座標 [m] real, intent(in) :: z(:) ! 鉛直座標 [m] real, intent(in) :: coril(size(r),size(z)) ! コリオリパラメータ [/s] real, intent(in) :: v(size(r),size(z)) ! 軸対称流 [m/s] real, intent(in) :: pres_s(size(z)) ! サウンディングの気圧 [Pa] real, intent(in) :: rho_s(size(z)) ! サウンディングの密度 [kg/m^3] real, intent(in), optional :: error ! イタレーションの収束条件 ! default = 1.0e-5 real, intent(inout) :: pres(size(r),size(z)) ! 平衡場の気圧 [Pa] real, intent(inout) :: rho(size(r),size(z)) ! 平衡場の密度 [kg/m^3] integer, intent(in), optional :: dl ! デバッグレベル real :: old_pres(size(r),size(z)), old_rho(size(r),size(z)) integer :: nr, nz integer :: i, j real :: err, err_tmp, err_max nr=size(r) nz=size(z) if(present(error))then err_max=error else err_max=1.0e-5 end if !-- 以下で各高度において, 密度は一定であるとして傾度風平衡から気圧を計算, !-- その値を用いて静力学平衡から密度を修正. eps 以下になるまで繰り返す. !-- 外縁で 2 次元場とサウンディングを一致. do i=1,nz old_pres(nr,i)=pres_s(i) end do !-- 密度については, 水平面一様で設定 do j=1,nz do i=1,nr old_rho(i,j)=rho_s(j) end do end do !-- 以下でイタレーション開始. err=err_max do while(err>=err_max) err=0.0 !-- 傾度風平衡から圧力場を計算 !$omp parallel default(shared) !$omp do schedule(dynamic) private(j) do j=1,nz call grad_wind_pres( r, coril(:,j), v(:,j), old_rho(:,j), r(nr), old_pres(nr,j), pres(:,j) ) end do !$omp end do !$omp end parallel !-- 静力学平衡から密度場を修正 do i=1,nr call grad_1d( z, pres(i,:), rho(i,:) ) do j=1,nz rho(i,j)=-rho(i,j)/g ! 静力学の式から, dp/dz=-g*rho であるので end do end do !-- 密度場の収束を計算 do j=1,nz do i=1,nr if(rho(i,j)==0.0)then err_tmp=abs(old_rho(i,j)-rho(i,j))/abs(old_rho(i,j)) else err_tmp=abs(old_rho(i,j)-rho(i,j))/abs(old_rho(i,j)) end if !-- 最大誤差の更新 if(err<=err_tmp)then err=err_tmp end if old_rho(i,j)=rho(i,j) old_pres(i,j)=pres(i,j) end do end do end do if(present(dl))then do j=1,nz call debug_flag_r( dl, 'typhoon_analy', 'hydro_grad_eqb (pres)', pres(1,j), 'Pa' ) end do end if end subroutine
Subroutine : | |||
x(2) : | real, intent(in)
| ||
y(2) : | real, intent(in)
| ||
c(2) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
nr : | integer, intent(inout)
| ||
stdopt : | logical, intent(in), optional
|
接線平均可能な半径を計算するルーチン
subroutine search_region_1d( x, y, c, r, nr, stdopt ) ! 接線平均可能な半径を計算するルーチン use Statistics implicit none real, intent(in) :: x(2) ! x 方向の両端座標 [m] real, intent(in) :: y(2) ! y 方向の両端座標 [m] real, intent(in) :: c(2) ! 円筒中心の位置座標 (x,y) [m] real, intent(in) :: r(:) ! 動径方向の位置座標 [m] integer, intent(inout) :: nr ! 平均可能半径の要素番号 (r(nr) が可能半径) logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) integer :: nrr, tmp_nr real :: xc, yc logical :: stderr nrr=size(r) xc=c(1) yc=c(2) nr=nrr if(present(stdopt))then stderr=stdopt else stderr=.false. end if if(abs(x(1)-xc) < r(nrr))then if(stderr.eqv..false.)then write(*,*) "typhoon_analy WARNING :" write(*,*) "|x(1)-xc| >= rmax. " write(*,*) "undef value is substituted out of region." end if call interpo_search_1d( r, abs(x(1)-xc), tmp_nr, stdopt=stderr ) nr=tmp_nr+1 ! interpo_search は abs の値より小さい r の要素番号が入るため. ! 以下も同じ理由 else if(abs(x(2)-xc) < r(nrr))then if(stderr.eqv..false.)then write(*,*) "typhoon_analy WARNING :" write(*,*) "|x(nx)-xc| >= rmax. " write(*,*) "undef value is substituted out of region." end if call interpo_search_1d( r, abs(x(2)-xc), tmp_nr, stdopt=stderr ) if(tmp_nr+1<nr)then nr=tmp_nr+1 end if else if(abs(y(1)-yc) < r(nrr))then if(stderr.eqv..false.)then write(*,*) "typhoon_analy WARNING :" write(*,*) "|y(1)-yc| >= rmax. " write(*,*) "undef value is substituted out of region." end if call interpo_search_1d( r, abs(y(1)-yc), tmp_nr, stdopt=stderr ) if(tmp_nr+1<nr)then nr=tmp_nr+1 end if else if(abs(y(2)-yc) < r(nrr))then if(stderr.eqv..false.)then write(*,*) "typhoon_analy WARNING :" write(*,*) "|y(ny)-yc| >= rmax. " write(*,*) "undef value is substituted out of region." end if call interpo_search_1d( r, abs(y(2)-yc), tmp_nr, stdopt=stderr ) if(tmp_nr+1<nr)then nr=tmp_nr+1 end if end if end if end if end if end subroutine search_region_1d
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
u(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
v(size(r),size(theta)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
| ||
axis : | character(2), intent(in), optional
|
任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. これ以外の場合, r で与えられる距離が領域外の範囲については undef で設定された値が代入される. undef がないときはゼロが代入される. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
から, 重線形内挿 interpolation_2d で計算.
本ルーチンは平面極座標グリッド値での偏差計算を行う. 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, tangent_mean_anom_scal_r2c を使用. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい.
subroutine tangent_conv_scal( x, y, xc, yc, u, r, theta, v, undef, undefg, undefgc, stdopt, axis ) ! 任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. ! これ以外の場合, r で与えられる距離が領域外の範囲については ! undef で設定された値が代入される. undef がないときはゼロが代入される. ! 平均化の手順は以下のとおり. ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値 ! から, 重線形内挿 interpolation_2d で計算. ! 本ルーチンは平面極座標グリッド値での偏差計算を行う. ! 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, ! tangent_mean_anom_scal_r2c を使用. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. use algebra use Derivation use max_min use statistics use Geometry use Map_Function implicit none real, intent(in) :: x(:) ! 右手座標系での第一成分 [m or rad] real, intent(in) :: y(:) ! 右手座標系での第二成分 [m or rad] real, intent(in) :: u(size(x),size(y)) ! 右手座標系での平均化する値 real, intent(in) :: xc ! 接線平均する際の中心 x 成分. [m or rad] real, intent(in) :: yc ! 接線平均する際の中心 y 成分. [m or rad] real, intent(in) :: r(:) ! (xc, yc) を中心とした極座標系動径座標 [m]. real, intent(in) :: theta(:) ! (xc, yc) を中心とした極座標系同位角座標 [rad]. real, intent(inout) :: v(size(r),size(theta)) ! 平均化した u のアノマリー. real, intent(in), optional :: undef ! 領域外の設定値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. ! デフォルトは "inc". logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) character(2), intent(in), optional :: axis ! x, y の座標系 ! 'xy' = デカルト座標系 [m] ! 'll' = 球面緯度経度座標系 [rad] ! デフォルトは 'xy'. integer :: i, j, nx, ny, nr, nt, i_undef, i_undefg real :: r_undef, r_undefg real :: work(size(r),size(theta)) real :: point(size(r),size(theta),2) integer :: ip(size(r),size(theta),2) real :: tmpx(2), tmpy(2), tmpz(2,2), inter(2) character(1) :: undefgcflag character(2) :: ax_flag logical, dimension(size(r)) :: undefgc_check logical :: ucf, stderr nx=size(x) ny=size(y) nr=size(r) nt=size(theta) undefgc_check(:)=.true. if(present(undef))then i_undef=int(undef) r_undef=undef else i_undef=-999 r_undef=-999.0 end if if(present(undefg))then i_undefg=int(undefg) r_undefg=undefg else i_undefg=-999 r_undefg=-999.0 end if if(present(undefg))then if(present(undefgc))then undefgcflag=undefgc(1:1) else undefgcflag="i" end if end if if(present(stdopt))then stderr=stdopt else stderr=.false. end if if(present(axis))then ax_flag=axis(1:2) else ax_flag='xy' end if !-- 先の引数条件をクリアしているか確認 --- !-- この操作は行わない. ! call search_region_1d( (/x(1),x(nx)/), (/y(1),y(ny)/), (/xc,yc/), r, nr, & ! & stdopt=stderr, axis=ax_flag(1:2) ) !-- 先に v に undef 値を入れておく. ! do j=1,nt ! do i=1,nrr ! v(i,j)=r_undef ! end do ! end do !-- 過程(1) --- if(ax_flag(1:2)=='xy')then do j=1,nt do i=1,nr call rt_2_xy( r(i), theta(j), point(i,j,1), point(i,j,2) ) point(i,j,1)=xc+point(i,j,1) point(i,j,2)=yc+point(i,j,2) end do end do else do j=1,nt do i=1,nr if(r(i)/=0.0)then call rt2ll( r(i), theta(j), xc, yc, point(i,j,1), point(i,j,2) ) else point(i,j,1)=xc point(i,j,2)=yc end if end do end do end if !-- 過程(2) --- do j=1,nt do i=1,nr call interpo_search_2d( x, y, point(i,j,1), point(i,j,2), ip(i,j,1), ip(i,j,2), undeff=i_undef, stdopt=stderr ) end do end do !-- 過程(3) --- do j=1,nt do i=1,nr if(ip(i,j,1)/=i_undef.and.ip(i,j,2)/=i_undef.and. ip(i,j,1)/=nx.and.ip(i,j,2)/=ny)then tmpx(1)=x(ip(i,j,1)) tmpx(2)=x(ip(i,j,1)+1) tmpy(1)=y(ip(i,j,2)) tmpy(2)=y(ip(i,j,2)+1) tmpz(1,1)=u(ip(i,j,1),ip(i,j,2)) tmpz(2,1)=u(ip(i,j,1)+1,ip(i,j,2)) tmpz(1,2)=u(ip(i,j,1),ip(i,j,2)+1) tmpz(2,2)=u(ip(i,j,1)+1,ip(i,j,2)+1) inter(1)=point(i,j,1) inter(2)=point(i,j,2) if(present(undefg))then ucf=undef_checker_2d( tmpz, undefg ) if(ucf.eqv..false.)then call interpolation_2d( tmpx, tmpy, tmpz, inter, work(i,j) ) else work(i,j)=r_undefg undefgc_check(i)=.false. end if else call interpolation_2d( tmpx, tmpy, tmpz, inter, work(i,j) ) end if else work(i,j)=r_undef end if end do end do do j=1,nt do i=1,nr v(i,j)=work(i,j) end do end do end subroutine tangent_conv_scal
Subroutine : | |||
signal : | character(2)
| ||
r(:) : | real, intent(in)
| ||
z(:) : | real, intent(in)
| ||
u(size(r),size(z)) : | real, intent(in)
| ||
v(size(r),size(z)) : | real, intent(in)
| ||
w(size(r),size(z)) : | real, intent(in)
| ||
rho(size(z)) : | real, intent(in)
| ||
nuh(size(r),size(z)) : | real, intent(in)
| ||
nuv(size(r),size(z)) : | real, intent(in)
| ||
val(size(r),size(z)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional | ||
sfctau(size(r)) : | real, intent(in), optional
|
円筒座標系におけるレイノルズ応力テンソルを計算する.
subroutine tangent_mean_Reynolds( signal, r, z, u, v, w, rho, nuh, nuv, val, undef, sfctau ) ! 円筒座標系におけるレイノルズ応力テンソルを計算する. use algebra use Derivation implicit none character(2) :: signal ! 計算するテンソル成分. ! ['11', '22', '33'] = それぞれ対角テンソル成分 ! ['12', '13', '21', '23', '31', '32'] = それぞれ非対角 ! テンソル成分. ただし, 対称テンソルであるため, '12'='21' を ! 計算していることに注意. real, intent(in) :: r(:) ! radial 方向の空間座標 [m] real, intent(in) :: z(:) ! vertical 方向の空間座標 [m] real, intent(in) :: u(size(r),size(z)) ! radial に対応する方向の 3 次元ベクトル成分 real, intent(in) :: v(size(r),size(z)) ! tangential に対応する方向の 3 次元ベクトル成分 real, intent(in) :: w(size(r),size(z)) ! vertical に対応する方向の 3 次元ベクトル成分 real, intent(in) :: rho(size(z)) ! 水平面に平均した基本場の密度 [kg/m^3] real, intent(in) :: nuh(size(r),size(z)) ! 水平渦粘性係数 real, intent(in) :: nuv(size(r),size(z)) ! 鉛直渦粘性係数 real, intent(inout) :: val(size(r),size(z)) ! 計算されたテンソル成分 ! 現在, 以下のオプションは使用していない. real, intent(in), optional :: undef real, intent(in), optional :: sfctau(size(r)) ! 地表面からのフラックス ! これが与えられれば, 最下層の応力はこれで置き換える. integer :: i ! イタレーション用添字 ! integer :: j ! イタレーション用添字 integer :: k ! イタレーション用添字 integer :: nr ! 空間配列要素数 1 次元目 integer :: nz ! 空間配列要素数 3 次元目 real :: dr ! 1 次元目を微分する格子間隔 [m] real :: dz ! 3 次元目を微分する格子間隔 [m] real :: sxx(size(r),size(z)), nu(size(r),size(z)) real :: stau(size(r)) dr=r(2)-r(1) dz=z(2)-z(1) nr=size(r) nz=size(z) val=0.0 stau=0.0 if(present(sfctau))then if(signal(2:2)=='3'.and.signal(1:1)/='3')then stau(:)=sfctau(:) end if end if !-- [NOTE] !-- 以下, 文字で case の or ができないため, !-- if 文の入れ子ではなく, if 文の並列表記で case と同じように見せかける. !-- これはもちろん, 上から順に if をたどるが, どの場合も 2 種類以上の if に !-- 合致しないことが既知であるために可能となる書き方であり, !-- 並列表記した if の 2 パターン以上に合致してしまうような条件文では, !-- case の代用には用いることができないことに注意. !-- 本ライブラリでこのような紛らわしい表記をしている場合は必ず NOTE が入る. if(signal(1:2)=='12'.or.signal(1:2)=='21')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if if(signal(1:2)=='23'.or.signal(1:2)=='32')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) if(signal(2:2)=='3')then do k=1,nz do i=1,nr nu(i,k)=nuv(i,k) end do end do else do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if end if if(signal(1:2)=='13'.or.signal(1:2)=='31')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) if(signal(2:2)=='3')then do k=1,nz do i=1,nr nu(i,k)=nuv(i,k) end do end do else do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if end if if(signal(1:2)=='11')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) call div( r, z, u, w, val ) do k=1,nz do i=1,nr if(r(i)/=0.0)then val(i,k)=val(i,k)+u(i,k)/r(i) end if end do end do do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if if(signal(1:2)=='22')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) call div( r, z, u, w, val ) do k=1,nz do i=1,nr if(r(i)/=0.0)then val(i,k)=val(i,k)+u(i,k)/r(i) end if end do end do do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if if(signal(1:2)=='33')then call tangent_mean_deform( signal, r, z, u, v, w, sxx ) call div( r, z, u, w, val ) do k=1,nz do i=1,nr if(r(i)/=0.0)then val(i,k)=val(i,k)+u(i,k)/r(i) end if end do end do do k=1,nz do i=1,nr nu(i,k)=nuh(i,k) end do end do end if !-- 以下の式は, 最初 val = 0 で if 文の中で計算されているものとされていない !-- ものに分かれるので, 統一の式で評価ができる. !-- 計算されていないものについてはそもそもゼロである. !-- 以下, 最下層は地表面フラックスを代入するかどうかのオプションのため, 別ループ if(present(sfctau))then do i=1,nr val(i,1)=stau(i) end do else do i=1,nr val(i,1)=rho(1)*nu(i,1)*(sxx(i,1)-(2.0/3.0)*val(i,1)) end do end if do k=2,nz do i=1,nr val(i,k)=rho(k)*nu(i,k)*(sxx(i,k)-(2.0/3.0)*val(i,k)) end do end do end subroutine
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
u(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
v(size(r),size(theta)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
| ||
axis : | character(2), intent(in), optional
|
任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. これ以外の場合, r で与えられる距離が領域外の範囲については undef で設定された値が代入される. undef がないときはゼロが代入される. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
から, 重線形内挿 interpolation_2d で計算.
(4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に偏差計算 Anomaly_1d 使用. 以上で各 nr について偏差値が得られる. 本ルーチンは平面極座標グリッド値での偏差計算を行う. 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, tangent_mean_anom_scal_r2c を使用. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい.
subroutine tangent_mean_anom_scal( x, y, xc, yc, u, r, theta, v, undef, undefg, undefgc, stdopt, axis ) ! 任意の物理量を台風の中心から接線方向へ平均し, そのアノマリーを計算するルーチン ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. ! これ以外の場合, r で与えられる距離が領域外の範囲については ! undef で設定された値が代入される. undef がないときはゼロが代入される. ! 平均化の手順は以下のとおり. ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値 ! から, 重線形内挿 interpolation_2d で計算. ! (4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に偏差計算 Anomaly_1d 使用. ! 以上で各 nr について偏差値が得られる. ! 本ルーチンは平面極座標グリッド値での偏差計算を行う. ! 接線平均を行ったのち, デカルト座標の偏差へ落とし込むには, ! tangent_mean_anom_scal_r2c を使用. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. use algebra use Derivation use max_min use statistics use Geometry implicit none real, intent(in) :: x(:) ! 右手座標系での第一成分 [m or rad] real, intent(in) :: y(:) ! 右手座標系での第二成分 [m or rad] real, intent(in) :: u(size(x),size(y)) ! 右手座標系での平均化する値 real, intent(in) :: xc ! 接線平均する際の中心 x 成分. [m or rad] real, intent(in) :: yc ! 接線平均する際の中心 y 成分. [m or rad] real, intent(in) :: r(:) ! (xc, yc) を中心とした極座標系動径座標 [m]. real, intent(in) :: theta(:) ! (xc, yc) を中心とした極座標系同位角座標 [rad]. real, intent(inout) :: v(size(r),size(theta)) ! 平均化した u のアノマリー. real, intent(in), optional :: undef ! 領域外の設定値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. ! デフォルトは "inc". logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) character(2), intent(in), optional :: axis ! x, y の座標系 ! 'xy' = デカルト座標系 [m] ! 'll' = 球面緯度経度座標系 [rad] ! デフォルトは 'xy'. integer :: i, j, nx, ny, nr, nt, i_undef, i_undefg real :: r_undef, r_undefg real :: work(size(r),size(theta)) character(3) :: undefgcflag character(2) :: ax_flag logical, dimension(size(r)) :: undefgc_check logical :: stderr nx=size(x) ny=size(y) nr=size(r) nt=size(theta) undefgc_check(:)=.true. if(present(undef))then i_undef=int(undef) r_undef=undef else i_undef=-999 r_undef=-999.0 end if if(present(undefg))then i_undefg=int(undefg) r_undefg=undefg else i_undefg=-999 r_undefg=-999.0 end if if(present(undefg))then if(present(undefgc))then undefgcflag=undefgc(1:3) else undefgcflag="inc" end if end if if(present(stdopt))then stderr=stdopt else stderr=.false. end if if(present(axis))then ax_flag=axis(1:2) else ax_flag='xy' end if !-- 過程(1) - (3) --- call tangent_conv_scal( x, y, xc, yc, u, r, theta, work, undef=r_undef, undefg=r_undefg, undefgc=undefgcflag(1:3), stdopt=stderr, axis=ax_flag(1:2) ) !-- 過程(4) --- if(present(undefg))then if(undefgcflag(1:1)=='i')then do i=1,nr call Anomaly_1d( work(i,:), v(i,:), error=undefg ) end do else do i=1,nr if(undefgc_check(i).eqv..true.)then call Anomaly_1d( work(i,:), v(i,:) ) else do j=1,nt v(i,j)=undefg end do end if end do end if else do i=1,nr call Anomaly_1d( work(i,:), v(i,:) ) end do end if end subroutine tangent_mean_anom_scal
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
scal(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
scal_anom(size(x),size(y)) : | real, intent(inout)
| ||
undef : | real, optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
|
台風中心から接線アノマリを計算し, デカルト座標系に戻す. 接線平均ルーチンを用いて物理量を平均し, その 1 次元データを確保しておく. 同時に, デカルト系での物理量の円筒座標系での radial 位置を求める. この radial の位置における接線平均値を先の 1 次元データから内挿で求める. この求めた内挿値を元のデカルトデータから引くことでアノマリとする. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい.
subroutine tangent_mean_anom_scal_Cart( x, y, xc, yc, scal, r, theta, scal_anom, undef, undefg, undefgc, stdopt ) ! 台風中心から接線アノマリを計算し, デカルト座標系に戻す. ! 接線平均ルーチンを用いて物理量を平均し, その 1 次元データを確保しておく. ! 同時に, デカルト系での物理量の円筒座標系での radial 位置を求める. ! この radial の位置における接線平均値を先の 1 次元データから内挿で求める. ! この求めた内挿値を元のデカルトデータから引くことでアノマリとする. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. use statistics implicit none real, intent(in) :: x(:) ! デカルト座標系での x 座標 real, intent(in) :: y(:) ! デカルト座標系での y 座標 real, intent(in) :: scal(size(x),size(y)) ! デカルト座標系での平均化する値. real, intent(in) :: xc ! 接線平均する際の中心 x 成分. real, intent(in) :: yc ! 接線平均する際の中心 y 成分. real, intent(in) :: r(:) ! 平均化したときの動径方向の座標(xc からの値を入れる). real, intent(in) :: theta(:) ! 平均化するときの接線方向の座標 [rad]. real, intent(inout) :: scal_anom(size(x),size(y)) ! デカルト系でのアノマリ. real, optional :: undef ! 内挿値が見つからないときの未定義値. ! デフォルトでは dcl の未定義値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) integer :: j, k, nx, ny, nr, nt, itmpr real :: tmp(size(r)) real :: tmpr, tmp_anom, undeff character(3) :: undefgcflag logical :: stderr if(present(undef))then undeff=undef else undeff=999.0 end if if(present(stdopt))then stderr=stdopt else stderr=.false. end if nx=size(x) ny=size(y) nr=size(r) nt=size(theta) if(present(undefg))then if(present(undefgc))then undefgcflag=undefgc(1:3) else undefgcflag="inc" end if end if if(present(undefg))then call tangent_mean_scal( x, y, xc, yc, scal, r, theta, tmp, undef=undeff, undefg=undefg, undefgc=undefgcflag(1:3), stdopt=stderr ) else call tangent_mean_scal( x, y, xc, yc, scal, r, theta, tmp, undef=undeff, stdopt=stderr ) end if !-- 接線平均値を内挿し, その内挿値を引いてアノマリを求める. if(present(undefg))then if(undefgcflag(1:3)=='err')then do k=1,ny do j=1,nx tmpr=sqrt((x(j)-xc)**2+(y(k)-yc)**2) call interpo_search_1d( r, tmpr, itmpr, stdopt=stderr ) if(nr>itmpr)then if(tmp(itmpr)/=undefg.and.tmp(itmpr+1)/=undefg)then call interpolation_1d( (/r(itmpr), r(itmpr+1)/), (/tmp(itmpr), tmp(itmpr+1)/), tmpr, tmp_anom ) scal_anom(j,k)=scal(j,k)-tmp_anom else scal_anom(j,k)=undefg end if else scal_anom(j,k)=undeff end if end do end do else do k=1,ny do j=1,nx tmpr=sqrt((x(j)-xc)**2+(y(k)-yc)**2) call interpo_search_1d( r, tmpr, itmpr, stdopt=stderr ) if(nr>itmpr)then if(tmp(itmpr)/=undefg.and.tmp(itmpr+1)/=undefg.and. scal(j,k)/=undefg)then call interpolation_1d( (/r(itmpr), r(itmpr+1)/), (/tmp(itmpr), tmp(itmpr+1)/), tmpr, tmp_anom ) scal_anom(j,k)=scal(j,k)-tmp_anom else scal_anom(j,k)=undefg end if else scal_anom(j,k)=undeff end if end do end do end if else do k=1,ny do j=1,nx tmpr=sqrt((x(j)-xc)**2+(y(k)-yc)**2) call interpo_search_1d( r, tmpr, itmpr, stdopt=stderr ) if(nr>itmpr)then call interpolation_1d( (/r(itmpr), r(itmpr+1)/), (/tmp(itmpr), tmp(itmpr+1)/), tmpr, tmp_anom ) scal_anom(j,k)=scal(j,k)-tmp_anom else scal_anom(j,k)=undeff end if end do end do end if end subroutine tangent_mean_anom_scal_Cart
Subroutine : | |||
charc : | character(6), intent(in)
| ||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
u1(size(x),size(y)) : | real, intent(in)
| ||
u2(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
v(size(r),size(theta)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
|
任意の物理量を台風の中心から接線方向へ平均アノマリを計算するルーチン 接線風速平均用. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. これ以外の場合, r で与えられる距離が領域外の範囲については undef で設定された値が代入される. undef がないときはゼロが代入される. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
から, 重線形内挿 interpolation_2d で計算. これと同じループ内で, 2 つのベクトル成分の平均値が得られるので, これらを用いて vec_prod_2d によって中心からの位置ベクトルとの外積を計算 . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.
(4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい. なお, 本ルーチンは可読性のため, 上の過程を tangent_mean_scal ルーチンに 丸投げしている. 本ルーチン内で行われるのは, charc に合わせて, 各デカルト座標点における中心点に対する接線風, 動径風を計算するのみである.
subroutine tangent_mean_anom_vec( charc, x, y, xc, yc, u1, u2, r, theta, v, undef, undefg, undefgc, stdopt ) ! 任意の物理量を台風の中心から接線方向へ平均アノマリを計算するルーチン ! 接線風速平均用. ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. ! これ以外の場合, r で与えられる距離が領域外の範囲については ! undef で設定された値が代入される. undef がないときはゼロが代入される. ! 平均化の手順は以下のとおり. ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値 ! から, 重線形内挿 interpolation_2d で計算. ! これと同じループ内で, 2 つのベクトル成分の平均値が得られるので, ! これらを用いて vec_prod_2d によって中心からの位置ベクトルとの外積を計算 ! . 同じループ内で中心からの距離で割って v の接線成分のみ抽出. ! (4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. ! 以上で各 nr について平均値が得られる. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. ! なお, 本ルーチンは可読性のため, 上の過程を tangent_mean_scal ルーチンに ! 丸投げしている. 本ルーチン内で行われるのは, charc に合わせて, ! 各デカルト座標点における中心点に対する接線風, 動径風を計算するのみである. use algebra use Derivation use max_min use statistics use Geometry implicit none character(6), intent(in) :: charc ! 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分. real, intent(in) :: x(:) ! デカルト座標系での x 座標 real, intent(in) :: y(:) ! デカルト座標系での y 座標 real, intent(in) :: u1(size(x),size(y)) ! デカルト座標系での平均化する値 1 real, intent(in) :: u2(size(x),size(y)) ! デカルト座標系での平均化する値 2 real, intent(in) :: xc ! 接線平均する際の中心 x 成分. real, intent(in) :: yc ! 接線平均する際の中心 y 成分. real, intent(in) :: r(:) ! 平均化したときの動径方向の座標(xc からの値を入れる). real, intent(in) :: theta(:) ! 平均化するときの接線方向の座標 [rad]. real, intent(inout) :: v(size(r),size(theta)) ! アノマリの u の値. real, intent(in), optional :: undef ! 領域外の設定値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) integer :: i, j, nx, ny, z_count real, dimension(size(x),size(y)) :: posx, posy, abpos, vecz logical :: stderr nx=size(x) ny=size(y) z_count=0 if(present(stdopt))then stderr=stdopt else stderr=.false. end if !-- まず, charc に合わせて, 水平風ベクトルを動径・接線風に変換する. !-- 変換すると, その値はスカラー値として評価されるので, !-- tangent_mean_scal に投げる. !-- 中心点に対する各デカルト座標点の水平位置ベクトルを計算. do j=1,ny do i=1,nx posx(i,j)=x(i)-xc posy(i,j)=y(j)-yc end do end do select case (trim(charc)) case ('vector') if(present(undefg))then call vec_prod_2d( posx, posy, u1, u2, vecz, undeff=undefg ) else call vec_prod_2d( posx, posy, u1, u2, vecz ) end if case ('scalar') if(present(undefg))then call dot_prod_2d( posx, posy, u1, u2, vecz, undeff=undefg ) else call dot_prod_2d( posx, posy, u1, u2, vecz ) end if case default write(*,*) "error : bad character. select 'vector', or 'scalar'." stop end select call abst_2d( posx, posy, abpos ) if(present(undefg))then do j=1,ny do i=1,nx if(vecz(i,j)/=undefg.and.abpos(i,j)/=0.0)then vecz(i,j)=vecz(i,j)/abpos(i,j) else if(abpos(i,j)==0.0)then z_count=z_count+1 vecz(i,j)=0.0 end if end do end do else do j=1,ny do i=1,nx if(abpos(i,j)/=0.0)then vecz(i,j)=vecz(i,j)/abpos(i,j) else z_count=z_count+1 vecz(i,j)=0.0 end if end do end do end if if(present(undef))then if(present(undefg))then call tangent_mean_anom_scal( x, y, xc, yc, vecz, r, theta, v, undef=undef, undefg=undefg, undefgc=trim(undefgc), stdopt=stderr ) else call tangent_mean_anom_scal( x, y, xc, yc, vecz, r, theta, v, undef=undef, stdopt=stderr ) end if else if(present(undefg))then call tangent_mean_anom_scal( x, y, xc, yc, vecz, r, theta, v, undefg=undefg, undefgc=trim(undefgc), stdopt=stderr ) else call tangent_mean_anom_scal( x, y, xc, yc, vecz, r, theta, v, stdopt=stderr ) end if end subroutine tangent_mean_anom_vec
Subroutine : | |||
signal : | character(2)
| ||
r(:) : | real, intent(in)
| ||
z(:) : | real, intent(in)
| ||
u(size(r),size(z)) : | real, intent(in)
| ||
v(size(r),size(z)) : | real, intent(in)
| ||
w(size(r),size(z)) : | real, intent(in)
| ||
val(size(r),size(z)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional |
デカルト座標系における変形速度テンソルを計算する.
subroutine tangent_mean_deform( signal, r, z, u, v, w, val, undef ) ! デカルト座標系における変形速度テンソルを計算する. use algebra use Derivation implicit none character(2) :: signal ! 計算するテンソル成分. ! ['11', '22', '33'] = それぞれ対角テンソル成分 ! ['12', '13', '21', '23', '31', '32'] = それぞれ非対角 ! テンソル成分. ただし, 対称テンソルであるため, '12'='21' を ! 計算していることに注意. real, intent(in) :: r(:) ! radial 方向の空間座標 [m] real, intent(in) :: z(:) ! vertical 方向の空間座標 [m] real, intent(in) :: u(size(r),size(z)) ! radial に対応する方向の 3 次元ベクトル成分 real, intent(in) :: v(size(r),size(z)) ! tangential に対応する方向の 3 次元ベクトル成分 real, intent(in) :: w(size(r),size(z)) ! vertical に対応する方向の 3 次元ベクトル成分 real, intent(inout) :: val(size(r),size(z)) ! 計算されたテンソル成分 ! 現在, 以下のオプションは使用していない. real, intent(in), optional :: undef integer :: i ! イタレーション用添字 integer :: j ! イタレーション用添字 integer :: k ! イタレーション用添字 integer :: nr ! 空間配列要素数 1 次元目 integer :: nz ! 空間配列要素数 2 次元目 real :: dr ! 1 次元目を微分する格子間隔 [m] real :: dz ! 2 次元目を微分する格子間隔 [m] dr=r(2)-r(1) dz=z(2)-z(1) nr=size(r) nz=size(z) !-- [NOTE] !-- 以下, 文字で case の or ができないため, !-- if 文の入れ子ではなく, if 文の並列表記で case と同じように見せかける. !-- これはもちろん, 上から順に if をたどるが, どの場合も 2 種類以上の if に !-- 合致しないことが既知であるために可能となる書き方であり, !-- 並列表記した if の 2 パターン以上に合致してしまうような条件文では, !-- case の代用には用いることができないことに注意. !-- 本ライブラリでこのような紛らわしい表記をしている場合は必ず NOTE が入る. if(signal(1:2)=='12'.or.signal(1:2)=='21')then do k=1,nz call grad_1d( r, v(:,k), val(:,k) ) do i=1,nr if(r(i)/=0.0)then val(i,k)=val(i,k)-v(i,k)/r(i) end if end do end do end if if(signal(1:2)=='23'.or.signal(1:2)=='32')then !$omp parallel default(shared) !$omp do schedule(runtime) private(k) do k=1,nr call grad_1d( z, v(k,:), val(k,:) ) end do !$omp end do !$omp end parallel end if if(signal(1:2)=='13'.or.signal(1:2)=='31')then call div( r, z, w, u, val ) end if if(signal(1:2)=='11')then !$omp parallel default(shared) !$omp do schedule(runtime) private(k) do k=1,nz call grad_1d( r, u(:,k), val(:,k) ) val(:,k)=2.0*val(:,k) end do !$omp end do !$omp end parallel end if if(signal(1:2)=='22')then !$omp parallel default(shared) !$omp do schedule(dynamic) private(j,k) do k=1,nz do j=1,nr if(r(j)/=0.0)then val(j,k)=2.0*u(j,k)/r(j) else val(j,k)=0.0 end if end do end do !$omp end do !$omp end parallel end if if(signal(1:2)=='33')then !$omp parallel default(shared) !$omp do schedule(runtime) private(j) do j=1,nr call grad_1d( z, w(j,:), val(j,:) ) val(j,:)=2.0*val(j,:) end do !$omp end do !$omp end parallel end if end subroutine
Subroutine : | |||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
u(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
v(size(r)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
| ||
axis : | character(2), intent(in), optional
|
任意の物理量を台風の中心から接線方向へ平均するルーチン このルーチンは接線風速を平均する時には用いることはできない. 接線の平均を行う際には, 別のルーチン, tangent_mean_vec の使用が必要. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. これ以外の場合, r で与えられる距離が領域外の範囲については undef で設定された値が代入される. undef がないときはゼロが代入される. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
から, 重線形内挿 interpolation_2d で計算.
(4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい.
subroutine tangent_mean_scal( x, y, xc, yc, u, r, theta, v, undef, undefg, undefgc, stdopt, axis ) ! 任意の物理量を台風の中心から接線方向へ平均するルーチン ! このルーチンは接線風速を平均する時には用いることはできない. ! 接線の平均を行う際には, 別のルーチン, tangent_mean_vec の使用が必要. ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. ! これ以外の場合, r で与えられる距離が領域外の範囲については ! undef で設定された値が代入される. undef がないときはゼロが代入される. ! 平均化の手順は以下のとおり. ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値 ! から, 重線形内挿 interpolation_2d で計算. ! (4) nr x nt 個の内挿スカラー値が求まったら, nt 方向に平均計算 mean_1d 使用. ! 以上で各 nr について平均値が得られる. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. use algebra use Derivation use max_min use statistics use Geometry implicit none real, intent(in) :: x(:) ! 右手座標系での第一成分 [m or rad] real, intent(in) :: y(:) ! 右手座標系での第二成分 [m or rad] real, intent(in) :: u(size(x),size(y)) ! 右手座標系での平均化する値 real, intent(in) :: xc ! 接線平均する際の中心 x 成分. [m or rad] real, intent(in) :: yc ! 接線平均する際の中心 y 成分. [m or rad] real, intent(in) :: r(:) ! (xc, yc) を中心とした極座標系動径座標 [m]. real, intent(in) :: theta(:) ! (xc, yc) を中心とした極座標系同位角座標 [rad]. real, intent(inout) :: v(size(r)) ! 平均化した u の値. real, intent(in), optional :: undef ! 領域外の設定値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. ! デフォルトは "inc". logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) character(2), intent(in), optional :: axis ! x, y の座標系 ! 'xy' = デカルト座標系 [m] ! 'll' = 球面緯度経度座標系 [rad] ! デフォルトは 'xy'. integer :: i, nx, ny, nr, nt, i_undef, i_undefg real :: r_undef, r_undefg real :: work(size(r),size(theta)) character(3) :: undefgcflag character(2) :: ax_flag logical, dimension(size(r)) :: undefgc_check logical :: stderr nx=size(x) ny=size(y) nr=size(r) nt=size(theta) undefgc_check(:)=.true. if(present(undef))then i_undef=int(undef) r_undef=undef else i_undef=-999 r_undef=-999.0 end if if(present(undefg))then i_undefg=int(undefg) r_undefg=undefg else i_undefg=-999 r_undefg=-999.0 end if if(present(undefg))then if(present(undefgc))then undefgcflag=undefgc(1:3) else undefgcflag="inc" end if end if if(present(stdopt))then stderr=stdopt else stderr=.false. end if if(present(axis))then ax_flag=axis(1:2) else ax_flag='xy' end if !-- 過程(1) - (3) call tangent_conv_scal( x, y, xc, yc, u, r, theta, work, undef=r_undef, undefg=r_undefg, undefgc=undefgcflag(1:3), stdopt=stderr, axis=ax_flag(1:2) ) !-- 過程(4) --- if(present(undefg))then if(undefgcflag(1:1)=='i')then do i=1,nr call Mean_1d( work(i,:), v(i), error=undefg ) end do else do i=1,nr if(undefgc_check(i).eqv..true.)then call Mean_1d( work(i,:), v(i) ) else v(i)=undefg end if end do end if else do i=1,nr call Mean_1d( work(i,:), v(i) ) end do end if end subroutine tangent_mean_scal
Subroutine : | |||
signal : | character(1)
| ||
r(:) : | real, intent(in)
| ||
z(:) : | real, intent(in)
| ||
u(size(r),size(z)) : | real, intent(in)
| ||
v(size(r),size(z)) : | real, intent(in)
| ||
w(size(r),size(z)) : | real, intent(in)
| ||
rho(size(z)) : | real, intent(in)
| ||
nuh(size(r),size(z)) : | real, intent(in)
| ||
nuv(size(r),size(z)) : | real, intent(in)
| ||
val(size(r),size(z)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional | ||
sfctau(size(r)) : | real, intent(in), optional
|
接線平均した乱流フラックスを計算する. 接線平均しているので, tau_{*2} 成分 (\theta 微分の成分) は含まれない.
subroutine tangent_mean_turb( signal, r, z, u, v, w, rho, nuh, nuv, val, undef, sfctau ) ! 接線平均した乱流フラックスを計算する. ! 接線平均しているので, tau_{*2} 成分 (\theta 微分の成分) は含まれない. use algebra use Derivation implicit none character(1) :: signal ! 円筒座標系の何番目の乱流成分かを判定する. ! [1] = 円筒座標における radial 座標成分 (方程式 vr 成分) ! [2] = 円筒座標における tangential 座標成分 (方程式 vt 成分) ! [3] = 円筒座標における vertical 座標成分 (方程式 w 成分) real, intent(in) :: r(:) ! 動径方向の位置座標 [m] real, intent(in) :: z(:) ! 鉛直方向の位置座標 [m] real, intent(in) :: u(size(r),size(z)) ! x に対応する方向の 2 次元ベクトル成分 real, intent(in) :: v(size(r),size(z)) ! y に対応する方向の 2 次元ベクトル成分 real, intent(in) :: w(size(r),size(z)) ! y に対応する方向の 2 次元ベ>クトル成分 real, intent(in) :: rho(size(z)) ! 水平面に平均した基本場の密度 [kg/m^3] real, intent(in) :: nuh(size(r),size(z)) ! 水平渦粘性係数 real, intent(in) :: nuv(size(r),size(z)) ! 鉛直渦粘性係数 real, intent(inout) :: val(size(r),size(z)) ! 乱流フラックス real, intent(in), optional :: undef real, intent(in), optional :: sfctau(size(r)) ! 地表面からのフラックス ! これが与えられれば, 最下層の応力はこれで置き換える. integer :: i ! イタレーション用添字 ! integer :: j ! イタレーション用添字 integer :: k ! イタレーション用添字 integer :: id ! イタレーション用添字 integer :: nr ! 空間配列要素数 1 次元目 integer :: nz ! 空間配列要素数 2 次元目 real :: dr ! 1 次元目を微分する格子間隔 [m] real :: dz ! 2 次元目を微分する格子間隔 [m] character(1) :: signaltau(3) real, dimension(size(r),size(z),3) :: tau ! signal 方向に ! 作用する 1,2,3 面に垂直な応力 real, dimension(size(r),size(z)) :: tmp real, dimension(size(r)) :: stau signaltau=(/ '1', '2', '3' /) dr=r(2)-r(1) dz=z(2)-z(1) nr=size(r) nz=size(z) val=0.0 do id=1,3 if(id/=2)then ! tau_{*2} 成分はゼロなので, 計算しない. if(present(sfctau))then stau(:)=sfctau(:) call tangent_mean_Reynolds( signal//signaltau(id), r, z, u, v, w, rho, nuh, nuv, tau(:,:,id), sfctau=stau ) else call tangent_mean_Reynolds( signal//signaltau(id), r, z, u, v, w, rho, nuh, nuv, tau(:,:,id) ) end if end if end do !-- (signal, 1) 成分の計算 do k=1,nz call grad_1d( r, tau(:,k,1), tmp(:,k)) do i=1,nr if(r(i)/=0.0)then val(i,k)=tmp(i,k)+val(i,k)+tau(i,k,1)/r(i) else val(i,k)=tmp(i,k)+val(i,k) end if end do end do !-- (signal, 3) 成分の計算 do i=1,nr call grad_1d( z, tau(i,:,3), tmp(i,:)) do k=1,nz val(i,k)=tmp(i,k)+val(i,k) end do end do end subroutine
Subroutine : | |||
charc : | character(6), intent(in)
| ||
x(:) : | real, intent(in)
| ||
y(:) : | real, intent(in)
| ||
xc : | real, intent(in)
| ||
yc : | real, intent(in)
| ||
u1(size(x),size(y)) : | real, intent(in)
| ||
u2(size(x),size(y)) : | real, intent(in)
| ||
r(:) : | real, intent(in)
| ||
theta(:) : | real, intent(in)
| ||
v(size(r)) : | real, intent(inout)
| ||
undef : | real, intent(in), optional
| ||
undefg : | real, intent(in), optional
| ||
undefgc : | character(3), intent(in), optional
| ||
stdopt : | logical, intent(in), optional
|
任意の物理量を台風の中心から接線方向へ平均するルーチン 接線風速平均用. 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. これ以外の場合, r で与えられる距離が領域外の範囲については undef で設定された値が代入される. undef がないときはゼロが代入される. 平均化の手順は以下のとおり. (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値
から, 重線形内挿 interpolation_2d で計算. これと同じループ内で, 2 つのベクトル成分の平均値が得られるので, これらを用いて vec_prod_2d によって中心からの位置ベクトルとの外積を計算 . 同じループ内で中心からの距離で割って v の接線成分のみ抽出.
(4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. 以上で各 nr について平均値が得られる. 以下, 処理の都合で所々に present(undefg) が入っているが, 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else の箇所を参照されたい. なお, 本ルーチンは可読性のため, 上の過程を tangent_mean_scal ルーチンに 丸投げしている. 本ルーチン内で行われるのは, charc に合わせて, 各デカルト座標点における中心点に対する接線風, 動径風を計算するのみである.
subroutine tangent_mean_vec( charc, x, y, xc, yc, u1, u2, r, theta, v, undef, undefg, undefgc, stdopt ) ! 任意の物理量を台風の中心から接線方向へ平均するルーチン ! 接線風速平均用. ! 引数の制限として, |x(1)-xc|, |xc-x(nx)|, |y(1)-yc|, |yc-y(ny)| > r. ! これ以外の場合, r で与えられる距離が領域外の範囲については ! undef で設定された値が代入される. undef がないときはゼロが代入される. ! 平均化の手順は以下のとおり. ! (1) nr, nt のすべての点についてそれに対応する x, y 座標値を rt_2_xy で計算. ! (2) その地点を含む x,y グリッドの微小領域を interpo_search_2d で検索. ! (3) その地点を含む 4 点が出たら, その地点でのスカラー値を 4 隅のスカラー値 ! から, 重線形内挿 interpolation_2d で計算. ! これと同じループ内で, 2 つのベクトル成分の平均値が得られるので, ! これらを用いて vec_prod_2d によって中心からの位置ベクトルとの外積を計算 ! . 同じループ内で中心からの距離で割って v の接線成分のみ抽出. ! (4) nr x nt 個の内挿接線成分値が求まったら, nt 方向に平均計算 mean_1d 使用. ! 以上で各 nr について平均値が得られる. ! 以下, 処理の都合で所々に present(undefg) が入っているが, ! 純粋な処理には関係ないので, ソースを読む場合は, present(undefg) の else ! の箇所を参照されたい. ! なお, 本ルーチンは可読性のため, 上の過程を tangent_mean_scal ルーチンに ! 丸投げしている. 本ルーチン内で行われるのは, charc に合わせて, ! 各デカルト座標点における中心点に対する接線風, 動径風を計算するのみである. use algebra use Derivation use max_min use statistics use Geometry implicit none character(6), intent(in) :: charc ! 動径成分か接線成分かの判別, vector = 接線, scalar = 動径成分. real, intent(in) :: x(:) ! デカルト座標系での x 座標 real, intent(in) :: y(:) ! デカルト座標系での y 座標 real, intent(in) :: u1(size(x),size(y)) ! デカルト座標系での平均化する値 1 real, intent(in) :: u2(size(x),size(y)) ! デカルト座標系での平均化する値 2 real, intent(in) :: xc ! 接線平均する際の中心 x 成分. real, intent(in) :: yc ! 接線平均する際の中心 y 成分. real, intent(in) :: r(:) ! 平均化したときの動径方向の座標(xc からの値を入れる). real, intent(in) :: theta(:) ! 平均化するときの接線方向の座標 [rad]. real, intent(inout) :: v(size(r)) ! 平均化した u の値. real, intent(in), optional :: undef ! 領域外の設定値 real, intent(in), optional :: undefg ! 格子点に欠損がある場合の内挿未定義値 character(3), intent(in), optional :: undefgc ! undefg がある場合の処理 ! "inc" = その格子点を参照値として内挿する点のみ平均操作時に除外して計算. ! "err" = その格子点を参照値として内挿する点を平均操作時に含む場合, 平均値そのものを未定義として計算. この場合, 未定義値は undefg となる. logical, intent(in), optional :: stdopt ! 探索範囲が見つからない旨の標準出力を表示させないようにする. ! default では .false. (表示させる) integer :: i, j, nx, ny, z_count real, dimension(size(x),size(y)) :: posx, posy, abpos, vecz logical :: stderr nx=size(x) ny=size(y) z_count=0 if(present(stdopt))then stderr=stdopt else stderr=.false. end if !-- まず, charc に合わせて, 水平風ベクトルを動径・接線風に変換する. !-- 変換すると, その値はスカラー値として評価されるので, !-- tangent_mean_scal に投げる. !-- 中心点に対する各デカルト座標点の水平位置ベクトルを計算. do j=1,ny do i=1,nx posx(i,j)=x(i)-xc posy(i,j)=y(j)-yc end do end do select case (trim(charc)) case ('vector') if(present(undefg))then call vec_prod_2d( posx, posy, u1, u2, vecz, undeff=undefg ) else call vec_prod_2d( posx, posy, u1, u2, vecz ) end if case ('scalar') if(present(undefg))then call dot_prod_2d( posx, posy, u1, u2, vecz, undeff=undefg ) else call dot_prod_2d( posx, posy, u1, u2, vecz ) end if case default write(*,*) "error : bad character. select 'vector', or 'scalar'." stop end select call abst_2d( posx, posy, abpos ) if(present(undefg))then do j=1,ny do i=1,nx if(vecz(i,j)/=undefg.and.abpos(i,j)/=0.0)then vecz(i,j)=vecz(i,j)/abpos(i,j) else if(abpos(i,j)==0.0)then z_count=z_count+1 vecz(i,j)=0.0 end if end do end do else do j=1,ny do i=1,nx if(abpos(i,j)/=0.0)then vecz(i,j)=vecz(i,j)/abpos(i,j) else z_count=z_count+1 vecz(i,j)=0.0 end if end do end do end if if(present(undef))then if(present(undefg))then call tangent_mean_scal( x, y, xc, yc, vecz, r, theta, v, undef=undef, undefg=undefg, undefgc=trim(undefgc), stdopt=stderr ) else call tangent_mean_scal( x, y, xc, yc, vecz, r, theta, v, undef=undef, stdopt=stderr ) end if else if(present(undefg))then call tangent_mean_scal( x, y, xc, yc, vecz, r, theta, v, undefg=undefg, undefgc=trim(undefgc), stdopt=stderr ) else call tangent_mean_scal( x, y, xc, yc, vecz, r, theta, v, stdopt=stderr ) end if end subroutine tangent_mean_vec
Function : | |||
undef_checker_1d : | logical | ||
val : | real, dimension(:), intent(in)
| ||
undef : | real, intent(in)
|
任意配列 val について, すべての要素について undef という値が入っているか どうかをチェックする. 1 つでも undef が入っていれば .true. を返す.
logical function undef_checker_1d( val, undef ) ! 任意配列 val について, すべての要素について undef という値が入っているか ! どうかをチェックする. 1 つでも undef が入っていれば .true. を返す. implicit none real, dimension(:), intent(in) :: val ! チェックする配列 real, intent(in) :: undef ! チェックする変数値 integer :: i, nx logical :: checker nx=size(val) checker=.false. do i=1,nx if(val(i)==undef)then checker=.true. exit end if end do undef_checker_1d=checker return end function
Function : | |||
undef_checker_2d : | logical | ||
val : | real, dimension(:,:), intent(in)
| ||
undef : | real, intent(in)
|
任意配列 val について, すべての要素について undef という値が入っているか どうかをチェックする. 1 つでも undef が入っていれば .true. を返す.
logical function undef_checker_2d( val, undef ) ! 任意配列 val について, すべての要素について undef という値が入っているか ! どうかをチェックする. 1 つでも undef が入っていれば .true. を返す. implicit none real, dimension(:,:), intent(in) :: val ! チェックする配列 real, intent(in) :: undef ! チェックする変数値 integer :: i, nx logical :: checker nx=size(val,2) checker=.false. do i=1,nx checker=undef_checker_1d( val(:,i), undef ) if(checker.eqv..true.)then exit end if end do undef_checker_2d=checker return end function
Function : | |||
undef_checker_3d : | logical | ||
val : | real, dimension(:,:,:), intent(in)
| ||
undef : | real, intent(in)
|
任意配列 val について, すべての要素について undef という値が入っているか どうかをチェックする. 1 つでも undef が入っていれば .true. を返す.
logical function undef_checker_3d( val, undef ) ! 任意配列 val について, すべての要素について undef という値が入っているか ! どうかをチェックする. 1 つでも undef が入っていれば .true. を返す. implicit none real, dimension(:,:,:), intent(in) :: val ! チェックする配列 real, intent(in) :: undef ! チェックする変数値 integer :: i, nx logical :: checker nx=size(val,3) checker=.false. do i=1,nx checker=undef_checker_2d( val(:,:,i), undef ) if(checker.eqv..true.)then exit end if end do undef_checker_3d=checker return end function