Class TimeFilter
In: util/timefilter.f90

Methods

Included Modules

GridSet

Public Instance methods

Function :
AsselinFilter(DimXMin:DimXMax, DimZMin:DimZMax) :real(8)
VarA(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
VarN(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(in)
VarB(DimXMin:DimXMax, DimZMin:DimZMax) :real(8), intent(inout)

時間フィルター; Asselin のタイムフィルターを利用

[Source]

  function AsselinFilter(VarA, VarN, VarB)
    !
    ! 時間フィルター; Asselin のタイムフィルターを利用
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in)     :: VarA(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(in)     :: VarN(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8), intent(inout)  :: VarB(DimXMin:DimXMax, DimZMin:DimZMax)  
    real(8)                 :: AsselinFilter(DimXMin:DimXMax, DimZMin:DimZMax)

    !時間フィルタ
    AsselinFilter = VarN + tfil * ( VarB  - 2.0d0 * VarN + VarA ) 
    
  end function AsselinFilter
Function :
xza_AsselinFilter_xza_xza_xza(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8)
xza_VarA(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(in)
xza_VarN(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(in)
xza_VarB(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8), intent(inout)

時間フィルター; Asselin のタイムフィルターを利用

[Source]

  function xza_AsselinFilter_xza_xza_xza( xza_VarA, xza_VarN, xza_VarB )
    !
    ! 時間フィルター; Asselin のタイムフィルターを利用
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8), intent(in)     :: xza_VarA(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8), intent(in)     :: xza_VarN(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8), intent(inout)  :: xza_VarB(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
    real(8)                 :: xza_AsselinFilter_xza_xza_xza(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)

    !時間フィルタ
    xza_AsselinFilter_xza_xza_xza = xza_VarN + tfil * ( xza_VarB  - 2.0d0 * xza_VarN + xza_VarA ) 
    
  end function xza_AsselinFilter_xza_xza_xza

[Validate]