Path: | src/dynamics/advectdiffusez.f90 |
Last Update: | Sat Apr 23 00:01:52 JST 2005 |
Copyright (C) GFD Dennou Club, 2004. All rights reserved.
* Developer: SUGIYAMA Ko-ichiro (sugiyama@gfd-dennou.org) * Version: $Id: advectdiffusez.f90,v 1.2 2005/04/22 15:01:52 sugiyama Exp $ * Tag Name: $Name: $ * Change History:
Z 方向の運動方程式中に現れる以下の項を 長い時間ステップ DelTimeLong において評価する.
* 移流項 * 浮力項 * 乱流拡散項 * 数値粘性項
* 移流項の計算には ***_adv となる変数を利用する. * 乱流拡散項, 数値粘性項の計算には ***_dif となる変数を利用する.
fs_VelX_adv(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in)
| ||
sf_VelZ_adv(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) | ||
ss_PotTemp_adv(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) | ||
fs_VelX_dif(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) | ||
sf_VelZ_dif(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) | ||
ss_Km_dif(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(in) | ||
sf_AdvDiffZ(DimXMin:DimXMax, DimZMin:DimZMax) : | real(8), intent(out)
|
subroutine AdvectDiffuseZ( fs_VelX_adv, sf_VelZ_adv, ss_PotTemp_adv, fs_VelX_dif, sf_VelZ_dif, ss_Km_dif, sf_AdvDiffZ ) !=begin !== Dependency !=end !== 暗黙の型宣言禁止 implicit none !=begin !== Input real(8), intent(in) :: fs_VelX_adv(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(in) :: sf_VelZ_adv(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(in) :: fs_VelX_dif(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(in) :: sf_VelZ_dif(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(in) :: ss_Km_dif(DimXMin:DimXMax, DimZMin:DimZMax) real(8), intent(in) :: ss_PotTemp_adv(DimXMin:DimXMax, DimZMin:DimZMax) !== Output real(8), intent(out) :: sf_AdvDiffZ(DimXMin:DimXMax, DimZMin:DimZMax) !=end !=== Work real(8) :: sf_AdvVelZ(DimXMin:DimXMax, DimZMin:DimZMax) real(8) :: sf_Buoy(DimXMin:DimXMax, DimZMin:DimZMax) real(8) :: sf_TurbVelZ(DimXMin:DimXMax, DimZMin:DimZMax) real(8) :: sf_NumDiffVelZ(DimXMin:DimXMax, DimZMin:DimZMax) call BeginSub("AdvectDiffuseZ", fmt="%c", c1="calculation of extra force term in z-direction.") !=== 移流項 sf_AdvVelZ = sf_avr_ff( ff_avr_fs( fs_VelX_adv ) * ff_dx_sf( sf_VelZ_adv ) ) + sf_avr_ss( ss_avr_sf( sf_VelZ_adv ) * ss_dz_sf( sf_VelZ_adv ) ) !=== 浮力項 sf_Buoy = sf_avr_ss( ( Grav * ss_PotTemp_adv / ss_PotTempBasicZ ) ) !=== 乱流拡散項 (時刻は t - Δt で評価) sf_TurbVelZ = sf_dx_ff( ff_avr_ss( ss_Km_dif ) * ff_dz_fs( fs_VelX_dif ) ) + sf_dx_ff( ff_avr_ss( ss_Km_dif ) * ff_dx_sf( sf_VelZ_dif ) ) + 2.0d0 * sf_dz_ss( ss_Km_dif * ss_dz_sf( sf_VelZ_dif ) ) - 2.0d0 * sf_dz_ss( ( ss_Km_dif ** 2.0d0 ) ) / ( 3.0d0 * ( Cm ** 2.0d0 ) * ( DelXZ ** 2.0d0 ) ) !=== 数値粘性項 (時刻は t - Δt で評価) sf_NumDiffVelZ = NuH * ( sf_dx_ff( ff_dx_sf( sf_VelZ_dif) ) ) + NuV * ( sf_dz_ss( ss_dz_sf( sf_VelZ_dif) ) ) !=== 外力項として出力 sf_AdvDiffZ = - sf_AdvVelZ + sf_TurbVelZ + sf_NumDiffVelZ + sf_Buoy call EndSub("AdvectDiffuseZ") end subroutine AdvectDiffuseZ