easyclimate.field.typhoon.axisymmetric

Axisymmetric Analysis

The axisymmetric analysis was conducted by spherical orthodrome transformation (Ritchie 1987; Nakamura et al. 1997; Yamazaki 2011). In the transformed coordinates, the cyclonic centre is relocated to the North Pole. The coordinate transform is easily conducted in three-dimensional Cartesian coordinates.

Consider a longitude–latitude grid, where the cyclone relocated at the North Pole (the NP coordinates hereafter). The transformation between the original and NP coordinates is straightforward with rotations in Cartesian coordinates. A positional vector in the Cartesian coordinates \(\pmb { r } = ( x , y , z ) ^ { \mathrm { T } }\) is obtained from the longitude– latitude coordinates by

\[x=\cos\lambda\sin\theta\]
\[y=\sin\lambda\sin\theta\]
\[z=\cos\theta\]

Now, \(\pmb { r }\) in the NP coordinates is projected to the original coordinates by following two rotations. Assuming that the cyclonic centre is \((\lambda_c, \theta_c)\), where \(\lambda_c\) and \(\theta_c\) are the central longitude and colatitude, respectively:

1. Rotate \(-\theta_c\) around the \(y\) axis (negative because of a clockwise rotation in the x–z plane). 2. Rotate \(\lambda_c\) around the \(z\) axis.

The two rotations are described by the following rotation matrix:

\[\begin{split}A _ { 1 } = \left[ \begin{array} { c c c } { { \cos \theta _ { \mathrm { c } } } } & { { 0 } } & { { \sin \theta _ { \mathrm { c } } } } \\ { { 0 } } & { { 1 } } & { { 0 } } \\ { { - \sin \theta _ { \mathrm { c } } } } & { { 0 } } & { { \cos \theta _ { \mathrm { c } } } } \end{array} \right]\end{split}\]
\[\begin{split}A _ { 2 } = \left[ \begin{array} { c c c } { \cos \lambda _ { \mathrm { c } } } & { - \sin \lambda _ { \mathrm { c } } } & { 0 } \\ { \sin \lambda _ { \mathrm { c } } } & { \cos \lambda _ { \mathrm { c } } } & { 0 } \\ { 0 } & { 0 } & { 1 } \end{array}\right]\end{split}\]
\[\begin{split}A = A_2 A_1 = \begin{bmatrix} \cos \lambda_c \cos \theta_c & -\sin \lambda_c & \cos \lambda_c \sin \theta_c \\ \sin \lambda_c \cos \theta_c & \cos \lambda_c & \sin \lambda_c \sin \theta_c \\ -\sin \theta_c & 0 & \cos \theta_c \end{bmatrix}\end{split}\]

The points in the NP coordinates are transformed to the original coordinates and are obtained by \(A\pmb{ r }\). The longitude and latitude corresponding to \(A\pmb{ r }\) are obtained by

\[\begin{split}\lambda = { \left\{ \begin{array} { l l } { \text{arctan} \! { \frac { y } { x } } { \bmod { 2 \pi } } , } & { x \! \neq \! 0 } \\ { 0 , } & { x \! = \! 0 } \end{array} \right. }\end{split}\]
\[\theta = a \cos z.\]

The scalar field is interpolated at (\(\lambda, \theta\)). The zonally symmetric and asymmetric components in the transformed coordinates represent the axially symmetric and asymmetric components, respectively.

The vectors such as winds need additional procedures. First, the horizontal winds expressed in the Cartesian coordinates are

\[\dot { x } = - u \sin \lambda - v \cos \lambda \cos \theta\]
\[\dot { y } = u \cos \lambda - v \sin \lambda \cos \theta\]
\[\dot { z } = v \sin \theta\]

and each component of a vector in the Cartesian coordinates is interpolated as a scalar. Then, the coordinates are rotated with \(A ^ { \mathrm { T } } \dot { x }\), where \(A^T = A_1 A_2\) , to obtain winds in the NP coordinates. The winds in the longitude–latitude coordinates are obtained from those in the Cartesian coordinates by

\[u = \dot { y } \cos \lambda - \dot { x } \sin \lambda\]
\[v = \operatorname{sgn} { ( \dot { z } ) } \sqrt { { ( \dot { x } \cos { \lambda } + \dot { y } \sin { \lambda } ) } ^ { 2 } + \dot { z } ^ { 2 } }\]

The transformed zonal and meridional winds represent tangential and radial components, respectively. In this study, the grid spacing is uniform in longitude (\(\lambda\)) and in colatitude (\(\theta=\dfrac{\pi}{2}-\phi\)), and the meridional extent is \(10^\circ\) from the North Pole.

See also

Functions

cyclone_axisymmetric_analysis(, polar_lat, lon_dim, ...)

Performs axisymmetric analysis of a cyclone by transforming data into a polar coordinate system centered on the cyclone.

Module Contents

easyclimate.field.typhoon.axisymmetric.cyclone_axisymmetric_analysis(data_input: xarray.DataArray, cyclone_center_point: Tuple[float, float], polar_lon: numpy.ndarray = np.arange(0, 360, 2), polar_lat: numpy.ndarray = np.arange(80, 90.1, 1), lon_dim: str = 'lon', lat_dim: str = 'lat', vertical_dim: str = 'level', R: float = 6371.0087714) easyclimate.core.datanode.DataNode

Performs axisymmetric analysis of a cyclone by transforming data into a polar coordinate system centered on the cyclone.

This function converts input data to a polar coordinate system based on the cyclone center, interpolates the data onto a polar grid, and decomposes it into symmetric and asymmetric components. The symmetric component is the azimuthal mean, and the asymmetric component is the deviation from this mean.

Parameters

data_inputxarray.DataArray

Input data with latitude, longitude, and optionally vertical dimensions, i.e., (lat, lon) or (level, lat, lon).

cyclone_center_pointTuple[float, float]

Cyclone center as (longitude, latitude) in degrees.

polar_lonnumpy.ndarray, optional

Array of longitudinal angles in degrees for the polar grid, default is np.arange(0, 360, 2).

polar_latnumpy.ndarray, optional

Array of latitudinal angles in degrees for the polar grid, default is np.arange(80, 90.1, 1).

lon_dimstr, optional

Name of the longitude dimension, default is ‘lon’.

lat_dimstr, optional

Name of the latitude dimension, default is ‘lat’.

vertical_dimstr, optional

Name of the vertical dimension, default is ‘level’.

Rfloat, optional ( \(\mathrm{km}\) ).

Earth’s radius, default is 6371.0087714.

Returns

easyclimate.DataNode

A DataNode containing three xarray.DataArray objects:

  • rotated: Data interpolated onto the polar grid.

  • rotated_symmetric: Azimuthal mean of the rotated data.

  • rotated_asymmetric: Deviation from the azimuthal mean.

See also

Example

>>> import xarray as xr
>>> import numpy as np
>>> data = xr.DataArray(np.random.rand(37, 241, 241), dims=['level', 'lat', 'lon'],
...                     coords={'level': np.array([1, 2, 3, 5, 7, 10, 20, 30, 50, 70, 100, 125, 150, 175, 200, 225, 250, 300, 350, 400, 450, 500, 550, 600, 650, 700, 750, 775, 800, 825, 850, 875, 900, 925, 950, 975, 1000]),
...                             'lat': np.arange(60, 0-0.25, -0.25),
...                             'lon': np.arange(110, 170 + 0.25, 0.25)
... )
>>> result = cyclone_axisymmetric_analysis(data, (140.20, 19.77))
>>> print(result)
<easyclimate.DataNode 'root'>
root: /
├── rotated: <xarray.DataArray>
│   Dimensions:  (level: 37, y: 11, polar_lon: 180)
│   Coordinates:
│     * lat        (y: 11): float64
│     * level      (level: 37): int32
│     * lon        (polar_lon: 180): int64
│     * polar_lat  (y: 11): float64
│     * polar_lon  (polar_lon: 180): int64
│     * y          (y: 11): float64
├── rotated_asymmetric: <xarray.DataArray>
│   Dimensions:  (level: 37, y: 11, polar_lon: 180)
│   Coordinates:
│     * lat        (y: 11): float64
│     * level      (level: 37): int32
│     * lon        (polar_lon: 180): int64
│     * polar_lat  (y: 11): float64
│     * polar_lon  (polar_lon: 180): int64
│     * y          (y: 11): float64
└── rotated_symmetric: <xarray.DataArray>
    Dimensions:  (level: 37, y: 11)
    Coordinates:
    * lat        (y: 11): float64
    * level      (level: 37): int32
    * polar_lat  (y: 11): float64
    * y          (y: 11): float64