GBVTD_main Module

The main module of GBVTD based on Lee et al. (1999, MWR)


Uses


Functions

public function matrix_sum(aij, akj, undeflag)

Calculate product for a component in a matrix

Arguments

Type IntentOptional Attributes Name
double precision, intent(in) :: aij(:)

Matrix A1

double precision, intent(in) :: akj(size(aij))

Matrix A2

logical, intent(in), optional :: undeflag(size(aij))

Undefined flag at each sampling point

Return Value doubleprecision


Subroutines

public subroutine Retrieve_velocity_GBVTD(nasym, r, t, td, Vd, RadTC, VT, VR, VT0, VR0, VTSn, VTCn, undef)

Solve unknown variables and return wind velocity on R-T coordinates based on the GBVTD technique.
------------------------------------------------------
[relationship between r and rh] --
------------------------------------------------------
i-1 i i+1
...|-- --|-- --|... : r(1:size(r)) = velocity radii

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nasym

wave number for asymmetric tangential wind

double precision, intent(in) :: r(:)

radial coordinate on which Vd is defined [m]

double precision, intent(in) :: t(:)

azimuthal coordinate on which Vd is defined [rad]

double precision, intent(in) :: td(size(r),size(t))

radar azimuthal angle defined at Vd(r,t) [rad]

double precision, intent(inout) :: Vd(size(r),size(t))

Doppler velocity defined on r-t [m s-1]

double precision, intent(in) :: RadTC

Distance from radar to TC center [m]

double precision, intent(out) :: VT(size(r),size(t))

retrieved total tangential wind [m s-1]

double precision, intent(out) :: VR(size(r),size(t))

retrieved total radial wind [m s-1]

double precision, intent(out) :: VT0(size(r),size(t))

retrieved axisymmetric radial component of rotating wind [m s-1]

double precision, intent(out) :: VR0(size(r),size(t))

retrieved axisymmetric tangential component of divergent wind [m s-1]

double precision, intent(out) :: VTSn(nasym,size(r),size(t))

retrieved tangential component of rotating wind [m s-1]

double precision, intent(out) :: VTCn(nasym,size(r),size(t))

retrieved radial component of rotating wind [m s-1]

double precision, intent(in), optional :: undef

undefined value for Vd

public subroutine calc_psid(rtc, r, theta, thetad, psid)

Calculate the nonlinear angle psid at a centain radius
From the geometry,
--(1)
--(2)
From Eq. (2),
--(3)
From Eqs. (1) and (3),

Arguments

Type IntentOptional Attributes Name
double precision, intent(in) :: rtc

Distance from the radar to TC center [m]

double precision, intent(in) :: r

Radius [m]

double precision, intent(in) :: theta(:)

azimuthal angle [deg]

double precision, intent(in) :: thetad(size(theta))

Azimuthal angle for radar [rad]

double precision, intent(out) :: psid(size(theta))

Nonlinear angle [rad]

public subroutine check_max(ival, mx_val)

Check maximum value in each element of the matrix "a"

Arguments

Type IntentOptional Attributes Name
double precision, intent(in) :: ival(:)

Matrix A

double precision, intent(out) :: mx_val

Max value in A