| Class | Thermo_Advanced_Routine |
| In: |
thermo_advanced_routine.f90
|
基礎ルーチン, 関数集を use して複雑な熱力学関数を計算するモジュール
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| dx : | real, intent(in)
| ||
| dy : | real, intent(in)
| ||
| dz : | real, intent(in)
| ||
| pt(nx,ny,nz) : | real, intent(in)
| ||
| BV(nx,ny,nz) : | real, intent(inout)
| ||
| undeff : | real, intent(in), optional |
ブラントバイサラ振動数を計算する
subroutine Brunt_Freq( nx, ny, nz, dx, dy, dz, pt, BV, undeff )
! ブラントバイサラ振動数を計算する
use analy
use Phys_Const
implicit none
integer, intent(in) :: nx ! 第 1 要素数
integer, intent(in) :: ny ! 第 2 要素数
integer, intent(in) :: nz ! 第 3 要素数
real, intent(in) :: dx ! x 方向の格子間隔 [m]
real, intent(in) :: dy ! y 方向の格子間隔 [m]
real, intent(in) :: dz ! z 方向の格子間隔 [m]
real, intent(in) :: pt(nx,ny,nz) ! 温位 [K]
real, intent(inout) :: BV(nx,ny,nz) ! ブラントバイサラ振動数 [1/s]
real, intent(in), optional :: undeff
integer :: i, j, k
if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
call grad_1d( nz, dz, pt(i,j,:), BV(i,j,:), undeff )
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
call grad_1d( nz, dz, pt(i,j,:), BV(i,j,:) )
end do
end do
!$omp end do
!$omp end parallel
end if
do k=1,nz
do j=1,ny
do i=1,nx
if(present(undeff))then
if(BV(i,j,k)==undeff)then
BV(i,j,k)=undeff
else
BV(i,j,k)=(g/pt(i,j,k))*BV(i,j,k)
end if
else
BV(i,j,k)=(g/pt(i,j,k))*BV(i,j,k)
end if
end do
end do
end do
end subroutine
| Subroutine : | |||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| dx : | real, intent(in)
| ||
| dy : | real, intent(in)
| ||
| dz : | real, intent(in)
| ||
| u(nx,ny,nz) : | real, intent(in)
| ||
| v(nx,ny,nz) : | real, intent(in)
| ||
| w(nx,ny,nz) : | real, intent(in)
| ||
| rho(nx,ny,nz) : | real, intent(in)
| ||
| pt(nx,ny,nz) : | real, intent(in)
| ||
| pv(nx,ny,nz) : | real, intent(inout)
| ||
| undeff : | real, intent(in), optional |
エルテルのポテンシャル渦度を計算する
subroutine Ertel_PV( nx, ny, nz, dx, dy, dz, u, v, w, rho, pt, pv, undeff )
! エルテルのポテンシャル渦度を計算する
use Thermo_Function
use Thermo_Routine
use analy
implicit none
integer, intent(in) :: nx ! 第 1 要素数
integer, intent(in) :: ny ! 第 2 要素数
integer, intent(in) :: nz ! 第 3 要素数
real, intent(in) :: dx ! x 方向の格子間隔 [m]
real, intent(in) :: dy ! y 方向の格子間隔 [m]
real, intent(in) :: dz ! z 方向の格子間隔 [m]
real, intent(in) :: u(nx,ny,nz) ! 速度場の x 成分 [m/s]
real, intent(in) :: v(nx,ny,nz) ! 速度場の y 成分 [m/s]
real, intent(in) :: w(nx,ny,nz) ! 速度場の z 成分 [m/s]
real, intent(in) :: rho(nx,ny,nz) ! 密度場 [kg/m^3]
real, intent(in) :: pt(nx,ny,nz) ! 温位場 [K]
real, intent(inout) :: pv(nx,ny,nz) ! PV [Km^2/kgs]
real, intent(in), optional :: undeff
real :: tmp1(nx,ny,nz)
real :: tmp2(nx,ny,nz)
real :: tmp3(nx,ny,nz)
real :: tmp4(nx,ny,nz)
real :: tmp5(nx,ny,nz)
real :: tmp6(nx,ny,nz)
real :: tmp7(nx,ny,nz)
integer :: i, j, k
if(present(undeff))then
! 温位の空間勾配を計算.
call grad_vec_3d( nx, ny, nz, dx, dy, dz, pt, tmp1, tmp2, tmp3, undeff )
! 3 次元 rotation を計算.
call rotate( nx, ny, nz, dx, dy, dz, u, v, w, tmp4, tmp5, tmp6, undeff )
! omega と grad pt の内積を計算
call dot_prod( tmp4, tmp5, tmp6, tmp1, tmp2, tmp3, tmp7, undeff )
else
! 温位の空間勾配を計算.
call grad_vec_3d( nx, ny, nz, dx, dy, dz, pt, tmp1, tmp2, tmp3 )
! 3 次元 rotation を計算.
call rotate( nx, ny, nz, dx, dy, dz, u, v, w, tmp4, tmp5, tmp6 )
! omega と grad pt の内積を計算
call dot_prod( tmp4, tmp5, tmp6, tmp1, tmp2, tmp3, tmp7 )
end if
if(present(undeff))then
!$omp parallel default(shared)
!$omp do private(i,j,k)
do k=1,nz
do j=1,ny
do i=1,nx
if(tmp7(i,j,k)==undeff.or.rho(i,j,k)==undeff)then
pv(i,j,k)=undeff
else
pv(i,j,k)=tmp7(i,j,k)/rho(i,j,k)
end if
end do
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do private(i,j,k)
do k=1,nz
do j=1,ny
do i=1,nx
pv(i,j,k)=tmp7(i,j,k)/rho(i,j,k)
end do
end do
end do
!$omp end do
!$omp end parallel
end if
end subroutine
| Subroutine : | |||
| z : | real, intent(in)
| ||
| z0m : | real, intent(in), dimension(:,:)
| ||
| richard : | real, intent(in), dimension(size(z0m,1),size(z0m,2))
| ||
| Lo : | real, intent(inout), dimension(size(z0m,1),size(z0m,2))
|
Louis(1980) で提案されている大気の不安定度を考慮したバルク係数の補正係数を計算する関数
subroutine Louis_horizon( z, z0m, richard, Lo )
! Louis(1980) で提案されている大気の不安定度を考慮したバルク係数の補正係数を計算する関数
use Thermo_Advanced_Function
implicit none
real, intent(in) :: z ! cm を求める高度 [m]
real, intent(in), dimension(:,:) :: z0m ! モデルで計算される粗度高度 [m]
real, intent(in), dimension(size(z0m,1),size(z0m,2)) :: richard ! バルクリチャードソン数
real, intent(inout), dimension(size(z0m,1),size(z0m,2)) :: Lo ! 補正係数
real, parameter :: b=5.0, c=5.0
real :: cm_tmp, zratio
integer :: i, j, nx, ny
nx=size(z0m,1)
ny=size(z0m,2)
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
Lo(i,j)=Louis( z, z0m(i,j), richard(i,j) )
end do
end do
!$omp end do
!$omp end parallel
end subroutine
| Subroutine : | |||
| za : | real, intent(in)
| ||
| pta : | real, intent(in), dimension(:,:)
| ||
| ptg : | real, intent(in), dimension(size(pta,1),size(pta,2))
| ||
| va : | real, intent(in), dimension(size(pta,1),size(pta,2))
| ||
| qva : | real, intent(in), dimension(size(pta,1),size(pta,2))
| ||
| qvs : | real, intent(in), dimension(size(pta,1),size(pta,2))
| ||
| Ri : | real, intent(inout), dimension(size(pta,1),size(pta,2))
|
バルクリチャードソン数を計算するルーチン
subroutine Rich_horizon( za, pta, ptg, va, qva, qvs, Ri )
! バルクリチャードソン数を計算するルーチン
use Thermo_Advanced_Function
implicit none
real, intent(in) :: za ! リチャードソン数を計算する高度 [m]
real, intent(in), dimension(:,:) :: pta ! za での仮温位 [K]
real, intent(in), dimension(size(pta,1),size(pta,2)) :: ptg ! 地表面での温位 [K]
real, intent(in), dimension(size(pta,1),size(pta,2)) :: va ! 高度 za での水平風速の絶対値 [m/s]
real, intent(in), dimension(size(pta,1),size(pta,2)) :: qva ! za での混合比 [kg/kg]
real, intent(in), dimension(size(pta,1),size(pta,2)) :: qvs ! 地表面での飽和混合比 [kg/kg]
real, intent(inout), dimension(size(pta,1),size(pta,2)) :: Ri ! 求めるリチャードソン数
real, dimension(size(pta,1),size(pta,2)) :: ptvg, ptva, dpt
integer :: i, j, nx, ny
nx=size(pta,1)
ny=size(pta,2)
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
Ri(i,j)=Rich( za, pta(i,j), ptg(i,j), va(i,j), qva(i,j), qvs(i,j) )
end do
end do
!$omp end do
!$omp end parallel
end subroutine
| Subroutine : | |||
| z : | real, intent(in)
| ||
| z0m : | real, intent(in), dimension(:,:)
| ||
| coem : | real, intent(inout), dimension(size(z0m,1),size(z0m,2))
| ||
| richard : | real, intent(in), dimension(size(z0m,1),size(z0m,2)), optional
|
運動量に関するバルク係数を計算するルーチン
subroutine cm_horizon( z, z0m, coem, richard )
! 運動量に関するバルク係数を計算するルーチン
use Thermo_Advanced_Function
implicit none
real, intent(in) :: z ! cm を求める高度 [m]
real, intent(in), dimension(:,:) :: z0m ! モデルで計算される粗度高度 [m]
real, intent(inout), dimension(size(z0m,1),size(z0m,2)) :: coem ! バルク係数
real, intent(in), dimension(size(z0m,1),size(z0m,2)), optional :: richard ! Louis (1980) のスキームで計算する場合のバルクリチャードソン数
integer :: i, j, nx, ny
nx=size(z0m,1)
ny=size(z0m,2)
if(present(richard))then
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
coem(i,j)=cm( z, z0m(i,j), richard(i,j) )
end do
end do
!$omp end do
!$omp end parallel
else
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
coem(i,j)=cm( z, z0m(i,j), richard(i,j) )
end do
end do
!$omp end do
!$omp end parallel
end if
end subroutine
| Subroutine : | |||
| cmd : | real, intent(in), dimension(:,:)
| ||
| va : | real, intent(in), dimension(size(cmd,1),size(cmd,2))
| ||
| velst : | real, intent(inout), dimension(size(cmd,1),size(cmd,2))
|
摩擦速度 u_* を計算するルーチン
subroutine ust_horizon( cmd, va, velst )
! 摩擦速度 u_* を計算するルーチン
use Thermo_Advanced_Function
implicit none
real, intent(in), dimension(:,:) :: cmd ! 高度 za でのバルク係数
real, intent(in), dimension(size(cmd,1),size(cmd,2)) :: va ! 高度 za での水平風の絶対値 [m/s]
real, intent(inout), dimension(size(cmd,1),size(cmd,2)) :: velst ! 摩擦速度 [m/s]
integer :: i, j, nx, ny
!$omp parallel default(shared)
!$omp do private(i,j)
do j=1,ny
do i=1,nx
velst(i,j)=ust( cmd(i,j), va(i,j) )
end do
end do
!$omp end do
!$omp end parallel
end subroutine