regrid2 package

Submodules

regrid2.crossSection module

class CrossSectionRegridder(latIn, latOut, levIn, levOut, latTypeIn=None, latSizeIn=None, latTypeOut=None, latSizeOut=None)[source]

Bases: object

PURPOSE: To perform all the tasks required to regrid the input data into the ouput data in the

latitude-level plane for all times

PROCEDURE:
Step One:

Make an instance of class CrossSectionRegridder passing it input and output grid information

Step Two:

Pass the input data with some descriptive parameters and get the output data in return

rgrd(dataIn, missingValueIn, missingMatch, logYes='yes', positionIn=None, maskIn=None, missingValueOut=None)[source]

To perform all the tasks required to regrid the input data, dataIn, into the ouput data, dataout in the latitude-level plane.

Parameters
dataIndata to regrid
missingValueInthe missing data value to use in setting missing in the mask. It is required.
  • None – there is no missing data

  • A number – if the value to use in the search for possible missing data. The presence of missing data at a grid point leads to recording 0.0 in the mask.

missingMatchthe comparison scheme used in searching for missing data in dataIn using the value

passed in as missingValueIn.

  • None – used if None is the entry for missingValueIn

  • exact – used if missingValue is the exact value from the file

  • greater – the missing data value is equal to or greater than missingValueIn

  • less – the missing data value is equal to or less than missingValueIn

logYeschoose the level regrid as linear in log of level or linear in level.
  • Set to ‘yes’ for log. Anything else is linear in level.

positionIna tuple with the numerical position of the dimensions
  • in C or Python order specified in the sequence latitude, level and time. Latitude and level are required. If time is missing submit None in its slot in the tuple. Notice that the length of the tuple is always three.

    Explicitly, in terms of the shape of dataIn as returned by python’s shape function

    positionIn[0] contains the position of latitude in dataIn positionIn[1] contains the position of level in dataIn or None positionIn[2] contains the position of time in dataIn or None

    If the c order shape of 3D data is (number of times, number of levels, number of latitudes) submit (2, 1, 0).

    If the c order shape of 2D data is (number of times, number of latitudes) submit (1, None, 0).

    Send in None if the shape is a subset of (time, level, latitude) which is evaluated as follows:

    • 2D – code assumes (1,0,None)

    • 3D – code assumes (2,1,0)

maskInan array of 1.0 and 0.0 values where the 0.0 value is used to mask the input data.

This mask only works on the latitude grid. It is not possible to mask out a region in the level plane. The 0.0 value removes the data from correponding grid point. The user can supply the following choices:

  • None – an array of 1.0s is created followed by substituting 0.0s for grid points with missing data in the input data array, dataIn

  • array – an array of 1.0s or 0.0s which must be either 2D or the actual size of the input data, dataIn. This user supplied mask might be used to mask a latitude region. It is not required to account for missing data in the input data. The code uses missingValueIn and missingMatch to supply the 0.0s for grid points with missing data in the input data array, dataIn.

missingValueOutthe value for the missing data used in writing the output data.
  • If left at the default entry, None, the code uses missingValueIn

  • If present or as a last resort 1.0e20

Returns
dataOutthe regridded data
checkdimension(x, name)[source]
Purposedimension checks
  1. has a len method

  2. data type is float32

  3. monotonically increasing vectors

Parameters
xcoordinate vector
namecoordinate vector ID
Returns
x, xsize – dimension vector and its size
generic_wts_bnds(lat)[source]
get_latitude_wts_bnds(checklatpass)[source]

Compare the passed checklatpass with the correct geophysical ones calculated here. After finding a match call the function to get the bounds.

wts,bnds = get_latitude_wts_bnds(checklatpass)

Parameters
checklatpassis the grid to check
Returns
wts, bnds - tuple with weights and bounds
get_region_latitude_wts_bnds(latRegionpass, latType, latSize)[source]

Routine : get_region_latitude_wts_bnds

Purposecompare the passed latitudes, latRegion, with the global

ones calculated here and extract the wts and bounds for

the region

Usagewts,bnds = get_region_latitude_wts_bnds(latRegion, latType, latSize)

where latRegion is the regional grid to check

Returns
wts, bnds - tuple with weights and bounds
latitude_bounds(lat_bnds)[source]

Purpose:

set up the shape and bounds for use by maparea

Usage:

Returns
tuple ( bn,bs )
rmserror(data1, data2)[source]

Purpose : compute the rms error for two data sets having the same shape

Passed : the two data sets

Returns
rms error
section(latvals, levvals)[source]

Purpose : make the crossi section analytical test case

Passed : the grid coordinate vectors

Returns
xsectiona temerature like cross section
sectionmask(dataIn, positionIn, maskIn, missingValueIn, missingMatch)[source]

Purpose : construct the mask for the input data for use by rgdlength

Usage : amskin = mask(dataIn, positionIn, maskIn, missingValueIn, missingValueOut, flag2D)

Returns
amskin
sendmsg(msg, value1=None, value2=None)[source]

Purpose : send the same message to the screen

Passed : msg - the string

Returns
return

regrid2.error module

exception RegridError(args='Unspecified error from regrid package')[source]

Bases: Exception

regrid2.esmf module

class EsmfRegrid(srcField, dstField, srcFrac=None, dstFrac=None, srcMaskValues=None, dstMaskValues=None, regridMethod=<RegridMethod.BILINEAR: 0>, ignoreDegenerate=False, unMappedAction=<UnmappedAction.IGNORE: 1>)[source]

Bases: object

Regrid source grid data to destination grid data

Constuct regrid object

Parameters
srcField :

the source field object of type EsmfStructFields

dstField :

the destination field object of type EsmfStructField

srcMaskValues :

Value of masked cells in source

dstMaskValues :

Value of masked cells in destination

srcFrac :

Cell fractions on source grid (type EsmfStructField

dstFrac :

Cell fractions on destination grid (type EsmfStructField)

regridMethod :

ESMF.RegridMethod.{BILINEAR,CONSERVE,PATCH}

unMappedAction :

ESMF.UnmappedAction.{IGNORE,ERROR}

ignoreDegenerate :

Ignore degenerate cells when checking inputs

getDstAreaFractions(rootPe)[source]

Get the destination grid fraction areas as used by conservative interpolation

Parameters
rootPeNone is local areas are returned, otherwise provide rootPe and the data will be gathered
Returns
numpy array
getDstAreas(rootPe)[source]

Get the dst grid areas as used by conservative interpolation

Parameters
rootPeNone is local areas are returned, otherwise provide rootPe and the data will be gathered
Returns
numpy array or None if interpolation is not conservative
getSrcAreaFractions(rootPe)[source]

Get the source grid fraction areas as used by conservative interpolation

Parameters
rootPeNone is local areas are returned, otherwise provide rootPe and the data will be gathered
Returns
numpy array
getSrcAreas(rootPe)[source]

Get the src grid areas as used by conservative interpolation

Parameters
rootPeNone is local areas are returned, otherwise provide rootPe and

the data will be gathered

Returns
numpy array or None if interpolation is not conservative
class EsmfStructField(esmfGrid, name, datatype, staggerloc=<StaggerLoc.CENTER: 0>)[source]

Bases: object

Structured field.

Creator for structured ESMF Field

Parameters
esmfGrid

instance of an ESMF

name field

name (must be unique)

datatype

data type, one of ‘float64’, ‘float32’, ‘int64’, or ‘int32’ (or equivalent numpy dtype)

staggerloc

ESMF.StaggerLoc.CENTER ESMF.StaggerLoc.CORNER

getData(rootPe)[source]

Get field data as a numpy array

Parameters
rootPeif None then local data will be fetched, otherwise gather the

data on processor “rootPe” (all other procs will return None).

Returns
numpy array or None.
getPointer()[source]

Get field data as a flat array

Returns
flat array pointer.
setLocalData(data, staggerloc, globalIndexing=False)[source]

Set local field data

Parameters
datafull numpy array, this method will take care of setting

a the subset of the data that reside on the local processor

staggerlocstagger location of the data
globalIndexingif True array was allocated over global index space, array

was allocated over local index space (on this processor)

class EsmfStructGrid(shape, coordSys=<CoordSys.SPH_DEG: 1>, periodicity=0, staggerloc=<StaggerLoc.CENTER: 0>, hasBounds=False)[source]

Bases: object

Structured grid

Parameters
shapeTuple of cell sizes along each axis
coordSyscoordinate system

ESMF.CoordSys.CART Cartesian ESMF.CoordSys.SPH_DEG (default) Degrees ESMF.CoordSys.SPH_RAD Radians

periodicityDoes the grid have a periodic coordinate

0 No periodicity 1 Periodic in x (1st) axis 2 Periodic in x, y axes

staggerlocESMF stagger location. ESMF.StaggerLoc.XXXX

The stagger constants are listed at the top

hasBoundsIf the grid has bounds, Run AddCoords for the bounds
getCellAreas()[source]

Get Cell Areas

Returns
cell areas or None if setCellAreas was not called
getCoordShape(staggerloc)[source]

Get the local coordinate shape (may be different on each processor)

Parameters
staggerloc(e.g. ESMF.StaggerLoc.CENTER)
Returns
tuple
getCoords(dim, staggerloc)[source]

Return the coordinates for a dimension

Parameters
dimdesired dimension (zero based indexing)
staggerlocStagger location
getLoHiBounds(staggerloc)[source]

Get the local lo/hi index values for the coordinates (per processor) (hi is not inclusive, lo <= index < hi)

Parameters
staggerloc(e.g. ESMF.StaggerLoc.CENTER)
Returns
lo, hi lists.
getLocalSlab(staggerloc)[source]

Get the local slab (ellipsis). You can use this to grab the data local to this processor

Parameters
staggerloc(e.g. ESMF.StaggerLoc.CENTER)
Returns
tuple of slices.
getMask(staggerloc=<StaggerLoc.CENTER: 0>)[source]

Get mask array. In ESMF, the mask is applied to cells.

Returns
mask numpy array 1 is invalid by default
setCellAreas(areas)[source]

Set the cell areas

Parameters
areasnumpy array
setCoords(coords, staggerloc=<StaggerLoc.CENTER: 0>, globalIndexing=False)[source]

Populate the grid with staggered coordinates (e.g. corner or center).

Parameters
coordsThe curvilinear coordinates of the grid. List of numpy arrays.

Must exist on all procs.

staggerlocThe stagger location ESMF.StaggerLoc.CENTER (default)

ESMF.StaggerLoc.CORNER

globalIndexingif True array was allocated over global index space, otherwise array was

allocated over local index space on this processor. This is only relevant if rootPe is None

Notes

coord dims in cdms2 are ordered in y, x, but ESMF expects x, y, hence the dimensions are reversed here.

setMask(mask, staggerloc=<StaggerLoc.CENTER: 0>)[source]

Set mask array. In ESMF, the mask is applied to cells.

Returns
mask numpy array 1 is invalid by default

Notes

This array exists on all procs

class EsmfUnstructGrid(numTopoDims, numSpaceDims)[source]

Bases: object

Constructor

Parameters
numTopoDimsnumber of topological dimensions
numSpaceDimsnumber of space dimensions
setCells(cellIndices, cellTypes, connectivity, cellMask=None, cellAreas=None)[source]

Set Cell connectivity.

Parameters
cellIndicesany 0-based.

cellTypes : any one of ESMF_MESHELEMTYPE_{TRI,QUAD,TETRA,HEX}.

connectivityNode : any connectivity array, see below for node ordering.

cellMask : any cellAreas area (volume) of each cell.

Notes

             3                       4-------------3
             /\                      |             |
            /  \                     |             |
           /    \                    |             |
          /      \                   |             |
         /        \                  |             |
        /          \                 |             |
       1------------2                1-------------2

          3                               8---------------7
         /|\                             /|              /|
        / | \                           / |             / |
       /  |  \                         /  |            /  |
      /   |   \                       /   |           /   |
     /    |    \                     5---------------6    |
    4-----|-----2                    |    |          |    |
     \    |    /                     |    4----------|----3
      \   |   /                      |   /           |   /
       \  |  /                       |  /            |  /
        \ | /                        | /             | /
         \|/                         |/              |/
          1                          1---------------2

ESMF_MESHELEMTYPE_TETRA             ESMF.MESHELEMTYPE_HEX
setNodes(indices, coords, peOwners=None)[source]

Set the nodal coordinates

Parameters
indicesIds of the nodes (0-based)
coordsnodal coordinates
peOwnersprocessor ranks where the coordinates reside (0-based)
toVTK(filename)[source]

Write grid to VTK file format

Parameters
filenameVTK file name

regrid2.git module

regrid2.gsRegrid module

Regridding of curvilinear structured grids

C_DOUBLE_P

alias of regrid2.gsRegrid.LP_c_double

class Regrid(src_grid, dst_grid, src_bounds=None, mkCyclic=False, handleCut=False, verbose=False)[source]

Bases: object

Constructor

Parameters
src_grid

source grid, a list of [x, y, …] coordinates or a cdms2.grid.Transient

dst_grid

destination grid, a list of [x, y, …] coordinate

src_bounds

list of cell bounding coordinates (to be used when handling a cut in coordinates)

mkCyclic

Add a column to the right side of the grid to complete a cyclic grid

handleCut

Add a row to the top of grid to handle a cut for grids such as the tri-polar grid verbose print diagnostic messages

apply(src_data_in, dst_data, missingValue=None)[source]

Apply interpolation

Parameters
src_data

data on source grid

dst_data

data on destination grid

missingValue

value that should be set for points falling outside the src domain, pass None if these should not be touched.

computeWeights(nitermax=100, tolpos=0.01)[source]

Compute the the interpolation weights

Parameters
nitermax

max number of iterations

tolpos

max tolerance when locating destination positions in index space

find(pattern, path)[source]
getDstGrid()[source]

Return the destination grid

Returns
grid
getIndicesAndWeights(dst_indices)[source]

Get the indices and weights for a single target location

Parameters
dst_indices

index set on the target grid

_: None

Returns
[index sets on original grid, weights]
getNumDstPoints()[source]

Return the number of points on the destination grid

Returns
number of points
getNumValid()[source]

Return the number of valid destination points. Destination points falling outside the source domain, more gnerally, points which could not be located on the source grid, reduce the number of valid points.

Returns
number of points
getPeriodicities()[source]

Get the periodicity lengths of the coordinates

Returns
numpy array, values inf indicate no periodicity
getSrcGrid()[source]

Return the source grid

Returns
grid
setMask(inDataOrMask)[source]

Set mask array. The mask is defined for nodes

Parameters
inDataOrMask
cdms2 array or flat mask array,

0 - valid data 1 - invalid data

_: None

Note: this definition is compatible with the numpy masked arrays
Note: note see setValidMask for the opposite definition
Note: should be called before computing the weights
setValidMask(inMask)[source]

Set valid mask array for the grid

Parameters
inMask
flat numpy array of type numpy.int32 or a valid cdms2 variable with its mask set.

0 - invalid, 1 - valid data

_: None

Note: This must be invoked before computing the weights, the mask is a property of the grid (not the data).
catchError(status, lineno)[source]
checkForCoordCut(coords, dims)[source]

Look for a cut in a coordinate system (e.g. tri-polar grid) Assume latitude is next to last coordinate and longitude is last coordinate!!!

Parameters
coords

input coordinates

dims

input dimensions

Returns
True for cut found

False for no cut

getTensorProduct(axis, dim, dims)[source]

Convert an axis into a curvilinear coordinate by applying a tensor product

Parameters
axis

1D array of coordinates

dim

dimensional index of the above coordinate

dims

sizes of all coordinates

Returns
coordinate values obtained by tensor product
handleCoordsCut(coords, dims, bounds)[source]

Generate connectivity across a cut. e.g. from a tri-polar grid. Assume latitude is next to last coordinate and longitude is last coordinate!!!

Parameters
coords

input coordinates list of rank

dims

input dimensions

bounds

boundaries for each coordinate

Returns
extended coordinates such that there is an extra row containing connectivity information across the cut
makeCoordsCyclic(coords, dims)[source]

Make coordinates cyclic

Parameters
coords

input coordinates

dims

input dimensions

Returns
new, extended coordinates such that the longitudes cover the sphere and new dimensions
makeCurvilinear(coords)[source]

Turn a mixture of axes and curvilinear coordinates into full curvilinear coordinates

Parameters
coords

list of coordinates

_: None

Returns
new list of coordinates and associated dimensions
test()[source]
testHandleCut()[source]
testMakeCyclic()[source]
testMasking()[source]

regrid2.gs_horizontal module

regrid2.horizontal module

class Horizontal(ingrid, outgrid)[source]

Bases: object

class Regridder(ingrid, outgrid)[source]

Bases: regrid2.horizontal.Horizontal

extractBounds(bounds)[source]
input_mask(ain, type, mask, missing=None)[source]

set up the input mask including missing from ain

regrid2.mvESMFRegrid module

ESMF regridding class

class ESMFRegrid(srcGridshape, dstGridshape, dtype, regridMethod, staggerLoc, periodicity, coordSys, srcGridMask=None, hasSrcBounds=False, srcGridAreas=None, dstGridMask=None, hasDstBounds=False, dstGridAreas=None, ignoreDegenerate=False, **args)[source]

Bases: regrid2.mvGenericRegrid.GenericRegrid

Regrid class for ESMF

Constructor

Parameters
srcGridShape

tuple source grid shape

dstGridShape

tuple destination grid shape

dtype

a valid numpy data type for the src/dst data

regridMethod

‘linear’, ‘conserve’, or ‘patch’

staggerLoc

the staggering of the field, ‘center’ or ‘corner’

periodicity

0 (no periodicity), 1 (last coordinate is periodic, 2 (both coordinates are periodic)

coordSys

‘deg’, ‘cart’, or ‘rad’

hasSrcBounds

tuple source bounds shape

hasDstBounds

tuple destination bounds shape

ignoreDegenerate

Ignore degenerate celss when checking inputs

apply(srcData, dstData, rootPe, globalIndexing=False, **args)[source]

Regrid source to destination. When used in parallel, if the processor is not the root processor, the dstData returns None.

Source data mask:

  • If you provide srcDataMask in args the source grid will be masked and weights will be recomputed.

  • Subsequently, if you do not provide a srcDataMask the last weights will be used to regrid the source data array.

  • By default, only the data are masked, but not the grid.

Parameters
srcDataarray source data, shape should cover entire global index space
dstDataarray destination data, shape should cover entire global index space
rootPeif other than None, then data will be MPI gathered on the specified rootPe processor
globalIndexingif True array was allocated over global index space, otherwise array was

allocated over local index space on this processor. This is only relevant if rootPe is None

args
computeWeights(**args)[source]

Compute interpolation weights

Parameters
args(not used)
fillInDiagnosticData(diag, rootPe)[source]

Fill in diagnostic data

Parameters
diaga dictionary whose entries, if present, will be filled valid

entries are: ‘srcAreaFractions’, ‘dstAreaFractions’, srcAreas’, ‘dstAreas’

rootPeroot processor where data should be gathered (or None if local areas are to be returned)
getDstAreaFractions(rootPe)[source]

Get the destination grid area fractions

Parameters
rootPeroot processor where data should be gathered (or None if local areas are to be returned)
Returns
fractional areas or None (if non-conservative)
getDstAreas(rootPe)[source]

Get the destination grid cell areas

Parameters
rootPeroot processor where data should be gathered (or None if local areas are to be returned)
Returns
areas or None if non-conservative interpolation
getDstGrid()[source]

Get the destination grid on this processor

Returns
grid
getDstLocalShape(staggerLoc)[source]

Get the local destination coordinate/data shape (may be different on each processor)

Parameters
staggerLoc(e.g. ‘center’ or ‘corner’)
Returns
tuple
getDstLocalSlab(staggerLoc)[source]

Get the destination local slab (ellipsis). You can use this to grab the data local to this processor

Parameters
staggerLoc(e.g. ‘center’)
Returns
tuple of slices
getSrcAreaFractions(rootPe)[source]

Get the source grid area fractions

Parameters
rootPeroot processor where data should be gathered (or None if local areas are to be returned)
Returns
fractional areas or None (if non-conservative)
getSrcAreas(rootPe)[source]

Get the source grid cell areas

Parameters
rootPeroot processor where data should be gathered (or None if local areas are to be returned)
Returns
areas or None if non-conservative interpolation
getSrcLocalShape(staggerLoc)[source]

Get the local source coordinate/data shape (may be different on each processor)

Parameters
staggerLoc(e.g. ‘center’ or ‘corner’)
Returns
tuple
getSrcLocalSlab(staggerLoc)[source]

Get the destination local slab (ellipsis). You can use this to grab the data local to this processor

Parameters
staggerLoc(e.g. ‘center’):
Returns
tuple of slices
setCoords(srcGrid, dstGrid, srcGridMask=None, srcBounds=None, srcGridAreas=None, dstGridMask=None, dstBounds=None, dstGridAreas=None, globalIndexing=False, **args)[source]

Populator of grids, bounds and masks

Parameters
srcGridlist [[z], y, x] of source grid arrays
dstGridlist [[z], y, x] of dstination grid arrays
srcGridMasklist [[z], y, x] of arrays
srcBoundslist [[z], y, x] of arrays
srcGridAreaslist [[z], y, x] of arrays
dstGridMasklist [[z], y, x] of array
dstBoundslist [[z], y, x] of arrays
dstGridAreaslist [[z], y, x] of arrays
globalIndexingif True array was allocated over global index space,

otherwise array was allocated over local index space on this processor. This is only relevant if rootPe is None

regrid2.mvGenericRegrid module

Generic interface to multiple regrid classes. No dependence on cdms2 variables.

class GenericRegrid(srcGrid, dstGrid, dtype, regridMethod, regridTool, srcGridMask=None, srcBounds=None, srcGridAreas=None, dstGridMask=None, dstBounds=None, dstGridAreas=None, **args)[source]

Bases: object

Generic Regrid class.

Constructor

Parameters
srcGridlist of numpy arrays, source horizontal coordinates
dstGrid

list of numpy arrays, destination horizontal coordinate

dtype

numpy data type for src/dst data

regridMethod

linear (bi, tri,…) default or conservative

regridTool

currently either ‘libcf’ or ‘esmf’

srcGridMask

array of same shape as srcGrid

srcBounds

list of numpy arrays of same shape as srcGrid

srcGridAreas

array of same shape as srcGrid

dstGridMask

array of same shape as dstGrid

dstBounds

list of numpy arrays of same shape as dstGrid

dstGridAreas

array of same shape as dstGrid

**args
additional arguments to be passed to the specific tool

libcf’: mkCyclic={True, False}, handleCut={True,False} ‘esmf’: periodicity={0,1}, coordSys={‘deg’, ‘cart’}, …

apply(srcData, dstData, rootPe=None, missingValue=None, **args)[source]

Regrid source to destination

Parameters
srcDataarray (input)
dstDataarray (output)
rootPeif other than None, then results will be MPI gathered
missingValueif not None, then data mask will be interpolated

and data value set to missingValue when masked

computeWeights(**args)[source]

Compute Weights

fillInDiagnosticData(diag, rootPe=None)[source]

Fill in diagnostic data

Parameters
diaga dictionary whose entries, if present, will be filled entries are tool dependent
rootPeroot processor where data should be gathered (or None if local areas are to be returned)
getDstGrid()[source]

Return the destination grid, may be different from the dst grid provided to the constructor due to domain decomposition

Returns
local grid on this processor
guessPeriodicity(srcBounds)[source]

Guess if a src grid is periodic

Parameters
srcBoundsthe nodal src set of coordinates
Returns
1 if periodic, warp around, 0 otherwise

regrid2.mvLibCFRegrid module

LibCF regridding class

class LibCFRegrid(srcGrid, dstGrid, srcGridMask=None, srcBounds=None, **args)[source]

Bases: regrid2.mvGenericRegrid.GenericRegrid

apply(srcData, dstData, missingValue=None, **args)[source]

Regrid source to destination

Parameters
srcData :

array (input)

dstData :

array (output)

missingValue :

value that should be set for points falling outside the src domain, pass None if these should not be touched.

computeWeights(**args)[source]

Compute interpolation weights

Parameters
**args arguments to be passed to gsRegrid, e.g. nitermax, tolpos, …
fillInDiagnosticData(diag, rootPe)[source]

Fill in diagnostic data

Parameters
diag

a dictionary whose entries, if present, will be filled valid entries are: ‘numDstPoints’ and ‘numValid’

rootPe

not used

getDstGrid()[source]

Get the grid of the dst data

Returns
grid
getSrcGrid()[source]

Get the grid of the src data (maybe larger than the dst grid passed to the constructor due to column/row padding)

Returns
grid

regrid2.mytest module

regrid2.pressure module

class PressureRegridder(axisIn, axisOut)[source]

Bases: object

PURPOSE: To perform all the tasks required to regrid the input data into the ouput data along

the pressure dimension only.

PROCEDURE:
Step One:

Make an instance of class PressureRegridder passing it input and output grid information

Step Two:

Pass the input data with some descriptive parameters and get the output data in return

rgrd(dataIn, missingValueIn, missingMatch, logYes='yes', positionIn=None, missingValueOut=None)[source]

To perform all the tasks required to regrid the input data, dataIn, into the ouput data, dataout along the level dimension only.

Parameters
dataIndata to regrid
missingValueInthe missing data value to use in setting missing in the mask. It is required.
  • None – there is no missing data

  • A number – if the value to use in the search for possible missing data. The presence of missing data at a grid point leads to recording 0.0 in the mask.

missingMatchthe comparison scheme used in searching for missing data in dataIn using the value
passed in as missingValueIn.
  • None – used if None is the entry for missingValueIn

  • exact – used if missingValue is the exact value from the file

  • greater – the missing data value is equal to or greater than missingValueIn

  • less – the missing data value is equal to or less than missingValueIn

logYeschoose the level regrid as linear in log of level or linear in level.
  • Set to ‘yes’ for log. Anything else is linear in level.

positionIna tuple with the numerical position of the dimensions
  • in C or Python order specified in the sequence longitude, latitude, level and time. Longitude and Latitude are required. If time is missing submit None in its slot in the tuple. Notice that the length of the tuple is always four.

    Explicitly, in terms of the shape of dataIn as returned by python’s shape function

    positionIn[0] contains the position of longitude in dataIn positionIn[1] contains the position of latitude in dataIn positionIn[2] contains the position of level in dataIn or None positionIn[3] contains the position of time in dataIn or None

    If the c order shape of 4D data is
    • (number of longitudes, number of times, number of levels, number of latitudes)

    submit
    • (0, 3, 2, 1).

    If the c order shape of 3D data is
    • (number of longitudes, number of times, number of latitudes)

    submit
    • (0, 2, 1, None).

    Send in None if the shape is a subset of (time, level, latitude, longitude) which is evaluated as follows:

    • 3D – code assumes (2,1,0,None)

    • 4D – code assumes (3,2,1,0)

missingValueOutthe value for the missing data used in writing the output data.
  • If left at the default entry, None, the code uses missingValueIn

  • If present or as a last resort 1.0e20

Returns
dataOutthe regridded data
checkorder(positionIn)[source]

Purpose :

construct the tuples for transposing the data to standard dimension order and the inverse for transposing it back to the original dimension order

Usage :

newOrder, inverseOrder = checkorder(positionIn)

Passed :

positionIn – array with location of longitude, latitude. level and time respectively in the sense of the python shape of the data

Returns
newOrdertuple to transpose data to the order (t,z,y,x)
inverseOrdertuple to transpose data to back to the original order
sendmsg(msg, value1=None, value2=None)[source]

Purpose :

send the same message to the screen

Passed :

msg - the string

value - the number associated with the string

Returns
return

regrid2.scrip module

class BicubicRegridder(outputGrid, remapMatrix, sourceAddress, destAddress, inputGrid=None, sourceFrac=None, destFrac=None)[source]

Bases: regrid2.scrip.ScripRegridder

Bicubic regrid.

Parameters
gradLat:

df/di

gradLon:

df/dj

gradLatlon:

d(df)/(di)(dj)

class BilinearRegridder(outputGrid, remapMatrix, sourceAddress, destAddress, inputGrid=None, sourceFrac=None, destFrac=None)[source]

Bases: regrid2.scrip.ScripRegridder

regrid(input)[source]
class ConservativeRegridder(outputGrid, remapMatrix, sourceAddress, destAddress, inputGrid=None, sourceFrac=None, destFrac=None, normalize='fracarea', normal=None, sourceArea=None, destArea=None)[source]

Bases: regrid2.scrip.ScripRegridder

First-order conservative regrid.

By default, the normalize option =”fracarea”, and array ‘normal’ is not specified.

If ‘normal’ is specified, it should be a one-dimensional array of the same length as the output grid size, with values:

1.0 for normalize=”fracarea”,

grid_frac for normalize=”destarea”, or

grid_frac*grid_area for normalize=”none”.

Parameters
sourceArea

is the area of the source grid cells

destArea

is the area of the destination grid cells

getDestinationArea()[source]
getSourceArea()[source]
regrid(input)[source]

call regridder

class DistwgtRegridder(outputGrid, remapMatrix, sourceAddress, destAddress, inputGrid=None, sourceFrac=None, destFrac=None)[source]

Bases: regrid2.scrip.ScripRegridder

regrid(input)[source]
class ScripRegridder(outputGrid, remapMatrix, sourceAddress, destAddress, inputGrid=None, sourceFrac=None, destFrac=None)[source]

Bases: object

getDestinationFraction()[source]
getInputGrid()[source]
getOutputGrid()[source]
getSourceFraction()[source]
readRegridder(fileobj, mapMethod=None, checkGrid=1)[source]

Read a regridder from an open fileobj.

Parameters
mapMethodis one of “conservative”, “bilinear”, “bicubic”, or “distwgt”.

If unspecified, it defaults to the method defined in the file.

If ‘checkGrid’ is 1 (default), the grid cells are checked for convexity, and ‘repaired’ if necessary.

Module contents

Interface to regridding facilities