
ilapsG
Inverts the Laplacian (on a gaussian grid) using spherical harmonics.
Prototype
function ilapsG ( zlap : numeric, zlmbda : numeric ) return_val [dimsizes(zlap)] : float or double
Arguments
zlapAn array of two or more dimensions containing the Laplacian values to invert. Last two dimensions must be lat and lon. Values must be in ascending latitudinal order. Grid must be global.
zlmbdaA constant if zlap is two-dimensional. If zlap has three or more dimensions then zlmbda must be the same size as zlap minus the rightmost two dimensions.
Return value
An array of the same size as zlap. Double if zlap or zlmbda is double, float otherwise.
Description
Inverts the Laplacian (on a gaussian grid) using spherical
harmonics. Missing values
are not allowed.
If zlmbda is identically zero, the poisson equation is solved. Otherwise,
the Helmholtz equation is solved.
Should not include the cyclic (wraparound) points. This is a restriction of the spherical harmonic (Spherepack) routines. Here is an example of passing in a two-dimensional array in which the last longitude point is left out.
z = sample ( x(:,0:nlon-2) ) ; does not include cyclic pointsIf zlap is on a fixed grid, ilapsF should be used. Also, note that ilapsG is the function version of ilapsg.
Use the ilapsG_Wrap function if metadata retention is desired. The interface is identical.
See Also
ilapsG_Wrap, ilapsg, ilapsF, ilapsf, lapsG, lapsF, ilapvf, ilapvg
Examples
Example 1
Read Z (on a gaussian grid) from a netCDF file and compute the inverse laplacian, solving the Poisson equation (k=0):
; These files are loaded by default in NCL V6.2.0 and newer load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_code.ncl" load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_csm.ncl" begin a = addfile("/cgd/cas/shea/MURPHYS/ATMOS/80.nc","r") z200 = a->Z(0,{189.},:,:) ; v200 is dimensioned nlat x nlon printVarSummary(z200) ilapl = ilapsG(z200,0) ; Use ilapsG_Wrap if metadata retention is desired ; ilapl = ilapsG_Wrap(z200,0) endExample 2
Read Z (on a gaussian grid) from a netCDF file and compute the inverse laplacian, solving the Poisson equation (k=0). Use ilapsG_Wrap to retain the metadata. Note that contributed.ncl needs to be loaded to use ilapsG_Wrap.
load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_code.ncl" load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/gsn_csm.ncl" load "$NCARG_ROOT/lib/ncarg/nclscripts/csm/contributed.ncl" begin a = addfile("/cgd/cas/shea/MURPHYS/ATMOS/80.nc","r") z200 = a->Z(0,{189.},:,:) ; z200 is dimensioned nlat x nlon printVarSummary(z200) ilapl = ilapsG(z200,0) ; Use ilapsG_Wrap if metadata retention is desired ; ilapl = ilapsG_Wrap(z200,0) end
Errors
If jer or ker is equal to:
1 : error in the specification of nlat
2 : error in the specification of nlon
4 : error in the specification of N (jer only)