Welcome to ESMPy - The ESMF Python Interface!
The ESMPy home page has all of the latest information on the ESMPy project including release notes, known bugs, supported platforms, and download information.
Please see the ESMF home page for more information on ESMF in general.
Fast Parallel Grid Remapping for Unstructured and Structured Grids gives a nice overview of the ESMF remapping functionality.
The ESMF Regridding Status page gives a good overview of the functionality that is available through various interfaces to ESMF regridding.
The ESMF_RegridWeightGen application is a command-line version of the functionality that is available through ESMPy.
Please contact esmf_support@list.woc.noaa.gov with any questions or problems.
ESMPy is a Python interface to the Earth System Modeling Framework (ESMF) regridding utility. ESMF is software for building and coupling weather, climate, and related models. ESMF has a robust, parallel and scalable remapping package, used to generate remapping weights. It can handle a wide variety of grids and options: logically rectangular grids and unstructured meshes; regional or global grids; 2D or 3D; and pole and masking options. ESMF also has capabilities to read grid information from NetCDF files in a variety of formats, including the evolving Climate and Forecast (CF) GridSpec and UGRID conventions.
ESMPy supports a single-tile logically rectangular discretization type called Grid and an unstructured discretization type called Mesh (ESMF also supports observational data streams). ESMPy supports bilinear, higher order patch recovery and first-order conservative regridding. There is also an option to ignore unmapped destination points and mask out points on either the source or destination. Regridding on the sphere takes place in 3D Cartesian space, so the pole problem is not an issue as it commonly is with many Earth system grid remapping softwares. Grid and Mesh objects can be created in 2D or 3D space, and 3D first-order conservative regridding is fully supported. Future plans for ESMPy involve the incorporation of observational data streams.
Regridding, also called remapping or interpolation, is the process of changing the grid underneath field data values while preserving the qualities of the original data. Different kinds of transformations are appropriate for different problems. Regridding may be needed when communicating data between Earth system model components such as land and atmosphere, or between different data sets to support analysis or visualization.
Regridding can be broken into two stages. The first stage is generation of an interpolation weight matrix that describes how points in the source grid contribute to points in the destination grid. The second stage is the multiplication of values on the source grid by the interpolation weight matrix to produce the appropriate values on the destination grid. ESMPy provides access to both stages through two separate interfaces.
There are many different interpolation methods, suitable for different problems. In ESMPy, the basic bilinear option is a two dimensional variant of linear interpolation. The higher order patch recovery is a second degree polynomial regridding method, which uses a least squares algorithm to calculate the polynomial. The first-order conservative regridding is a variant of a constant method which compares the proportions of overlapping source and destination cells to determine appropriate weights. All of these methods can be broken down to a simple sparse matrix multiplication operation between interpolation weights and data values.
The ESMF User’s Guide contains information on building and installing ESMF. The ESMF Reference Manual contains information on the architecture of ESMF, example code, and details of the API (Application Programming Interface).
Instructions on how to download the ESMPy code can be found at the ESMPy Download page.
The following packages are required to work with ESMPy:
The following packages are optional:
mpi4py - python bindings to MPI, needed to run the parallel regridding tests
nose - for nose testing
Installation of ESMPy requires a pointer to a file named esmf.mk inside of an ESMF installation. This file resides in a directory which looks like:
<ESMF_INSTALL_DIR>/lib/lib<g<or>O>/<platform>/esmf.mk
If the ESMFMKFILE flag is set when building ESMPy then it will not need to be referenced again. If not, an environment variable of the same name must be set with the path to the esmf.mk file EVERY time that a new shell is initiated.
The ESMPy build can be installed in a custom location using the –prefix, –home, or –install-base flags to the install command. If this is done, then this location needs to be added to the PYTHONPATH environment variable in the user’s shell EVERY time that a new shell is initiated. If a customized install location is not specified, ESMPy will be installed in the standard Python package installation directory on that particular machine.
Note: The ESMPy build does not have to be installed to be used. The PYTHONPATH environment variable can simply be pointed to the directory containing the ESMF module (esmfcontrib-ESMPy/src from a default git clone) after the build command.
As usual, any command followed by –help should print out some information on what options are available.
An installation of ESMPy in the default location for Python packages can be done with the following command issued from the top level ESMPy directory:
default Python package installation:
python setup.py build –ESMFMKFILE=<DIR_TO_esmf.mk> install
custom install location:
python setup.py build –ESMFMKFILE=<DIR_TO_esmf.mk>
python setup.py install –prefix=<custom_install_location>
setenv PYTHONPATH <custom_install_location>/lib/*/site_packages
Please contact esmf_support@list.woc.noaa.gov with any questions or problems.
The setup.py file can be used to run all of the ESMPy tests (in serial or parallel), like this:
python setup.py test_all
python setup.py test_all_parallel
or subsets of the tests individually, like this:
python setup.py test
python setup.py test_regrid
python setup.py test_regrid_from_file
python setup.py test_parallel
python setup.py test_regrid_parallel
python setup.py test_regrid_from_file_parallel
NOTE: The regrid_from_file tests can take up a lot of memory and bandwidth. The “test_regrid_from_file_dryrun” command will simply download the test files without actually running them (allowing the stress on the machine to be applied to bandwidth first, and then memory).
Alternatively, if the nose package is available, the tests can be run with:
nosetests
Individual tests can be run with nose using the following format:
nosetests <file>:<test>
e.g.
nosetests src/ESMF/test/unit_test.py:field_regrid_test
ESMPy doesn’t include many aspects of ESMF, including components, array interfaces, time management, etc. The limitations listed here are relative to ESMF offline and integrated regridding capabilities.
Testing related:
ESMPy can create Grid or Mesh objects from NetCDF files in a variety of formats. A Mesh can be created from files in SCRIP, ESMF, and UGRID formats. Grid files can be in SCRIP and GRIDSPEC format.
This file format is used by the SCRIP [4] package, grid files that work with that package should also work here. SCRIP format files are capable of storing either 2D logically rectangular grids or 2D unstructured grids. More information can be found in the ESMF reference manual section on the SCRIP Grid File Format.
ESMF has custom unstructured grid file format for describing meshes. This format is more compatible than the SCRIP format with the methods used to create a Mesh object, so less conversion needs to be done to create a Mesh. The ESMF format is thus more efficient than SCRIP when used with ESMPy. More information can be found in the ESMF reference manual section on the ESMF Unstructured Grid File Format.
GRIDSPEC is an extension to the Climate and Forecast (CF) metadata conventions for the representation of gridded data for Earth System Models. ESMPy supports NetCDF files that follow the CF GRIDSPEC convention to support logically rectangular lat/lon grids. More information can be found in the ESMF reference manual section on the CF Convention GRIDSPEC File Format.
UGRID is an extension to the CF metadata conventions for the unstructured grid data model. ESMPy support NetCDF files that follow the CF UGRID convention for unstructured grids. More information can be found in the ESMF reference manual section on the CF Convention UGRID File Format.
When creating a Mesh from a SCRIP format file, there are a number of options to control the output Mesh. The data is located at the center of the grid cell in a SCRIP grid; whereas the data is located at the corner of a cell in an ESMF Mesh object. Therefore, we create a Mesh object by default by constructing a “dual” mesh using the coordinates in the file. If the user wishes to not construct the dual mesh, the optional argument ‘convert_to_dual’ may be used to control this behavior. When ‘convert_to_dual’ is set to False, the Mesh constructed from the file will not be the dual. This is necessary when the Mesh is part of a conservative regridding operation, so the weights are properly generated for the cell centers in the file.
A Mesh may also be created with boolean flags to specify whether or not to add an area property to the Mesh ‘add_user_area’, or to add a mask ‘add_mask’ held by the NetCDF variable indicated in the optional argument, ‘varname’. These argument are only valid for UGRID formatted files.
A number of optional boolean arguments are also supported to create a structured Grid from a file. These include ‘is_sphere’ to indicate whether the grid is spherical or regional, ‘add_corner_stagger’ to add the corner stagger information to the Grid for conservative regridding, and ‘add_user_area’ to specify whether to read in the cell area from the NetCDF file or to calculate them. Also, for GRIDSPEC formmated files there is the ‘add_mask’ optional argument to add a mask held by the NetCDF variable indicated in optional argument, ‘varname’, and the ‘coord_names’ argument to specify the longitude and latitude variable names in GRIDSPEC file containing multiple sets of coordinates.
The following three sections describe the regridding methods that are available in ESMPy.
In 2D, ESMPy supports bilinear regridding between any combination of the following:
In 3D, ESMPy supports bilinear regridding between any combination of the following:
To use the bilinear method the user may created their Fields on any stagger location for Grids or the node location (MeshLoc.NODE) for Meshes. For Grids, the stagger location upon which the Field was built must contain coordinates.
In 2D, ESMPy supports patch regridding between any combination of the following:
Patch regridding is currently not supported in 3D.
To use the patch method the user may created their Fields on any stagger location for Grids or the node location (MeshLoc.NODE) for Meshes. For Grids, the stagger location upon which the Field was built must contain coordinates.
See references [1] and [2] for more information.
In 2D, ESMPy supports first-order conservative regridding between any combination of the following:
In 3D, ESMPy supports first-order conservative regridding between any combination of the following:
To use the first-order conservative method the user must have created their Fields on the center stagger location (StaggerLoc.CENTER in 2D or StaggerLoc.CENTER_VCENTER in 3D) for Grids or the element location (MeshLoc.ELEMENT) for Meshes. For Grids, the corner stagger location (StaggerLoc.CORNER in 2D or StaggerLoc.CORNER_VFACE in 3D) must contain coordinates describing the outer perimeter of the Grid cells.
See reference [3] for more information.
Masking is the process whereby parts of a grid can be marked to be ignored during an operation, such as regridding. Masking can be used on a source grid to indicate that certain portions of the grid should not be used to generate regridded data. This is useful, for example, if a portion of a source grid contains unusable values. Masking can also be used on a destination grid to indicate that the portion of the field built on that part of the grid should not receive regridded data. This is useful, for example, when part of the grid isn’t being used (e.g. the land portion of an ocean grid).
ESMPy currently supports masking for Fields built on structured Grids and element masking for Fields built on unstructured Meshes. A Grid mask is initialized by setting mask values in the Numpy Array returned from the Grid.get_item() call using the ‘item’ variable. A Mesh mask is initialized by passing mask values into the Mesh.add_elements() call using the ‘element_mask’ variable. The Field mask can then be setup by indicating the values to use for the mask in the ‘mask_values’ variable of the Field constructor. However, the Field mask does not need to be setup to mask values in the regridding operation. Regrid masking is handled by passing the mask values into the ‘src_mask_values’ or ‘dst_mask_values’ variables of the Regrid constructor. For example, if ‘dst_mask_values’ is set to (/1,2/), then any location in the Grid or Mesh that has a value of 1 or 2 will be masked.
Masking behavior differs slightly between regridding methods. For non-conservative regridding methods (e.g. bilinear or high-order patch), masking is done on points. For these methods, masking a destination point means that the point won’t participate in regridding (e.g. won’t receive an interpolated value). For these methods, masking a source point means that the entire source cell using that point is masked out. In other words, if any corner point making up a source cell is masked then the whole cell is masked. For conservative regridding methods (e.g. first-order conservative) masking is done on cells. Masking a destination cell means that the cell won’t participate in regridding (e.g. won’t receive an interpolated value). Similarly, masking a source cell means that the cell won’t participate in regridding (e.g. won’t contribute to interpolation). For any type of interpolation method (conservative or non-conservative) the masking is set on the location upon which the Fields passed into the regridding call are built. For example, if Fields built on StaggerLoc.CENTER are passed into the Regrid() call then the masking should also be set in StaggerLoc.CENTER.
In the case that the Grid is on a sphere (coord_sys=CoordSys.SPH_DEG or CoordSys.SPH_RAD) then the coordinates given in the Grid are interpreted as latitude and longitude values. The coordinates can either be in degrees or radians as indicated by the ‘coord_sys’ flag set during Grid creation. As is true with many global models, this application currently assumes the latitude and longitude refer to positions on a perfect sphere, as opposed to a more complex and accurate representation of the earth’s true shape such as would be used in a GIS system.
If a destination point cannot be mapped to a location in the source grid, the user has two options. The user may ignore those destination points that cannot be mapped by setting the ‘unmapped_action’ argument to UnmappedAction.IGNORE. The user also has the option to return an error if unmapped destination points exist. This is the default behavior, so the user can either not set the ‘unmapped_action’ argument or the user can set it to UnmappedAction.ERROR. At this point ESMPy does not support extrapolation to destination points outside the unmasked source Field.
Class | Description |
---|---|
Manager | A manager class to initialize and finalize ESMF |
Field | A data field built on a Grid or Mesh |
Grid | A structured grid for coordinate representation |
Mesh | An unstructured grid for coordinate representation |
Regrid | The regridding utility |
Named constants | Description |
---|---|
CoordSys | Specify the coordinate system of a Grid |
DecompFlag | Specify how DistGrid elements are decomposed over DEs |
FileFormat | Specify the format of a data file |
GridItem | Specify a mask or area item on a Grid |
LogKind | Specify how much logging should be done |
MeshElemType | Specify the type of the Mesh elements |
MeshLoc | Specify a nodal or elemental Mesh |
Region | Specify various regions in the data layout of |
RegridMethod | Specify which interpolation method to use regridding |
StaggerLoc | Specify the position for data in a Grid cell |
TypeKind | Specify the type and kind of data |
UnmappedAction | Specify which action to take with respect to unmapped destination points |
PoleMethod | Specify which type of artificial pole to construct on the source Grid for regridding |
The first few tutorials are stand-alone scripts that can be run from any Python interpeter.
import ESMF # Start up ESMF, this call is only necessary to override the default parameters # for logkind (ESMF.LogKind.NONE) and debug (False) esmpy = ESMF.Manager(logkind=ESMF.LogKind.MULTI, debug=True) print "Hello ESMPy World from PET (processor) {0}!".format(ESMF.local_pet())
# This example demonstrates how to create an ESMPy grid from file # The grid file is required, it can be retrieved from the ESMF data repository: # wget http://www.earthsystemmodeling.org/download/data/ll2.5deg_grid.nc import ESMF # Start up ESMF, this call is only necessary to override the default parameters # for logkind (ESMF.LogKind.NONE) and debug (False) esmpy = ESMF.Manager(logkind=ESMF.LogKind.MULTI, debug=True) # Create a uniform global latlon grid from a SCRIP formatted file grid = ESMF.Grid(filename="ll2.5deg_grid.nc", filetype=ESMF.FileFormat.SCRIP) # Create a field on the centers of the grid field = ESMF.Field(grid, "field", staggerloc=ESMF.StaggerLoc.CENTER) print "Successfully read a grid and created a field!" print "The field values on PET (processor) # {0} are:".format(ESMF.local_pet()) print field
# This example demonstrates how to create an ESMPy mesh from file # The mesh file is required, it can be retrieved from the ESMF data repository: # wget http://www.earthsystemmodeling.org/download/data/mpas_uniform_10242_dual_counterclockwise.nc import ESMF # Start up ESMF, this call is only necessary to override the default parameters # for logkind (ESMF.LogKind.NONE) and debug (False) esmpy = ESMF.Manager(logkind=ESMF.LogKind.MULTI, debug=True) # create an ESMF formatted unstructured mesh with clockwise cells removed mesh = ESMF.Mesh(filename="mpas_uniform_10242_dual_counterclockwise.nc", filetype=ESMF.FileFormat.ESMFMESH) # create a field on the nodes of the mesh field = ESMF.Field(mesh, "field", meshloc=ESMF.MeshLoc.NODE) print "Successfully read a mesh and created a field!" print "The field values on PET (processor) # {0} are:".format(ESMF.local_pet()) print field
# This example demonstrates how to regrid between a grid and a mesh. # The grid and mesh files are required, they can be retrieved from the ESMF data repository: # wget http://www.earthsystemmodeling.org/download/data/ll2.5deg_grid.nc # wget http://www.earthsystemmodeling.org/download/data/mpas_uniform_10242_dual_counterclockwise.nc import ESMF grid1 = "data/ll2.5deg_grid.nc" grid2 = "data/mpas_uniform_10242_dual_counterclockwise.nc" #SCRIP # Create a uniform global latlon grid from a SCRIP formatted file grid = ESMF.Grid(filename="ll2.5deg_grid.nc", filetype=ESMF.FileFormat.SCRIP, add_corner_stagger=True) # create a field on the center stagger locations of the source grid srcfield = ESMF.Field(grid, 'srcfield', staggerloc=ESMF.StaggerLoc.CENTER) # initialize the field to a constant value srcfield[...] = 25 # create an ESMF formatted unstructured mesh with clockwise cells removed mesh = ESMF.Mesh(filename="mpas_uniform_10242_dual_counterclockwise.nc", filetype=ESMF.FileFormat.ESMFMESH) # create a field on the elements of the destination mesh dstfield = ESMF.Field(mesh, 'dstmesh', meshloc=ESMF.MeshLoc.ELEMENT) # create an object to regrid data from the source to the destination field regrid = ESMF.Regrid(srcfield, dstfield, \ regrid_method=ESMF.RegridMethod.CONSERVE, \ unmapped_action=ESMF.UnmappedAction.IGNORE) # do the regridding from source to destination field dstfield = regrid(srcfield, dstfield) print "Successfully read a grid and a mesh and did a regridding of a constant field!" print "The field values on PET (processor) # {0} are:".format(ESMF.local_pet()) print dstfield.data
The following tutorials will show you how to build all of the pieces necessary to regrid data between Fields built on Grids and Meshes.
def grid_create(bounds, coords, domask=False, doarea=False): ''' PRECONDITIONS: 'bounds' contains the number of indices required for the two dimensions of a 2D Grid. 'coords' contains the upper and lower coordinate bounds of the Grid. 'domask' is a boolean value that gives the option to put a mask on this Grid. 'doarea' is an option to create user defined areas on this Grid. \n POSTCONDITIONS: A 2D Grid has been created. RETURN VALUES: \n Grid :: grid \n ''' lb_x = float(bounds[0]) lb_y = float(bounds[1]) ub_x = float(bounds[2]) ub_y = float(bounds[3]) min_x = float(coords[0]) min_y = float(coords[1]) max_x = float(coords[2]) max_y = float(coords[3]) cellwidth_x = (max_x-min_x)/(ub_x-lb_x) cellwidth_y = (max_y-min_y)/(ub_y-lb_y) cellcenter_x = cellwidth_x/2 cellcenter_y = cellwidth_y/2 max_index = np.array([ub_x,ub_y]) grid = ESMF.Grid(max_index, coord_sys=ESMF.CoordSys.CART) ## CORNERS grid.add_coords(staggerloc=[ESMF.StaggerLoc.CORNER]) # get the coordinate pointers and set the coordinates [x,y] = [0,1] gridCorner = grid.coords[ESMF.StaggerLoc.CORNER] for i in xrange(gridCorner[x].shape[x]): gridCorner[x][i, :] = float(i)*cellwidth_x + \ min_x + grid.lower_bounds[ESMF.StaggerLoc.CORNER][x] * cellwidth_x # last line is the pet specific starting point for this stagger and dim for j in xrange(gridCorner[y].shape[y]): gridCorner[y][:, j] = float(j)*cellwidth_y + \ min_y + grid.lower_bounds[ESMF.StaggerLoc.CORNER][y] * cellwidth_y # last line is the pet specific starting point for this stagger and dim ## CENTERS grid.add_coords(staggerloc=[ESMF.StaggerLoc.CENTER]) # get the coordinate pointers and set the coordinates [x,y] = [0,1] gridXCenter = grid.get_coords(x, staggerloc=ESMF.StaggerLoc.CENTER) gridYCenter = grid.get_coords(y, staggerloc=ESMF.StaggerLoc.CENTER) for i in xrange(gridXCenter.shape[x]): gridXCenter[i, :] = float(i)*cellwidth_x + cellwidth_x/2.0 + \ min_x + grid.lower_bounds[ESMF.StaggerLoc.CENTER][x] * cellwidth_x # last line is the pet specific starting point for this stagger and dim for j in xrange(gridYCenter.shape[y]): gridYCenter[:, j] = float(j)*cellwidth_y + cellwidth_y/2.0 + \ min_y + grid.lower_bounds[ESMF.StaggerLoc.CENTER][y] * cellwidth_y # last line is the pet specific starting point for this stagger and dim if domask: mask = grid.add_item(ESMF.GridItem.MASK) mask[:] = 1 mask[np.where((1.75 <= gridXCenter.data < 2.25) & (1.75 <= gridYCenter.data < 2.25))] = 0 if doarea: grid.add_item(ESMF.GridItem.AREA) area = grid.get_item(ESMF.GridItem.AREA) area[:] = 5.0 return grid
def grid_create_3d(bounds, coords, domask=False, doarea=False): ''' PRECONDITIONS: 'bounds' contains the number of indices required for the two dimensions of a 2D Grid. 'coords' contains the upper and lower coordinate bounds of the Grid. 'domask' is a boolean value that gives the option to put a mask on this Grid. 'doarea' is an option to create user defined areas on this Grid. \n POSTCONDITIONS: A Grid has been created. RETURN VALUES: \n Grid :: grid \n ''' lb_x = float(bounds[0]) lb_y = float(bounds[1]) lb_z = float(bounds[2]) ub_x = float(bounds[3]) ub_y = float(bounds[4]) ub_z = float(bounds[5]) min_x = float(coords[0]) min_y = float(coords[1]) min_z = float(coords[2]) max_x = float(coords[3]) max_y = float(coords[4]) max_z = float(coords[5]) cellwidth_x = (max_x-min_x)/(ub_x-lb_x) cellwidth_y = (max_y-min_y)/(ub_y-lb_y) cellwidth_z = (max_z-min_z)/(ub_z-lb_z) cellcenter_x = cellwidth_x/2 cellcenter_y = cellwidth_y/2 cellcenter_z = cellwidth_z/2 max_index = np.array([ub_x,ub_y,ub_z]) grid = ESMF.Grid(max_index, coord_sys=ESMF.CoordSys.CART) ## CORNERS grid.add_coords(staggerloc=[ESMF.StaggerLoc.CORNER_VFACE]) # get the coordinate pointers and set the coordinates [x,y,z] = [0,1,2] gridCorner = grid.coords[ESMF.StaggerLoc.CORNER_VFACE] for i in xrange(gridCorner[x].shape[x]): gridCorner[x][i, :, :] = float(i)*cellwidth_x + \ min_x + grid.lower_bounds[ESMF.StaggerLoc.CORNER_VFACE][x] * cellwidth_x # last line is the pet specific starting point for this stagger and dim for j in xrange(gridCorner[y].shape[y]): gridCorner[y][:, j, :] = float(j)*cellwidth_y + \ min_y + grid.lower_bounds[ESMF.StaggerLoc.CORNER_VFACE][y] * cellwidth_y # last line is the pet specific starting point for this stagger and dim for k in xrange(gridCorner[z].shape[z]): gridCorner[z][:, :, k] = float(k)*cellwidth_z + \ min_z + grid.lower_bounds[ESMF.StaggerLoc.CORNER_VFACE][z] * cellwidth_z # last line is the pet specific starting point for this stagger and dim ## CENTERS grid.add_coords(staggerloc=[ESMF.StaggerLoc.CENTER_VCENTER]) # get the coordinate pointers and set the coordinates [x,y,z] = [0,1,2] gridXCenter = grid.get_coords(x, staggerloc=ESMF.StaggerLoc.CENTER_VCENTER) gridYCenter = grid.get_coords(y, staggerloc=ESMF.StaggerLoc.CENTER_VCENTER) gridZCenter = grid.get_coords(z, staggerloc=ESMF.StaggerLoc.CENTER_VCENTER) for i in xrange(gridXCenter.shape[x]): gridXCenter[i, :, :] = float(i)*cellwidth_x + cellwidth_x/2.0 + \ min_x + grid.lower_bounds[ESMF.StaggerLoc.CENTER_VCENTER][x] * cellwidth_x # last line is the pet specific starting point for this stagger and dim for j in xrange(gridYCenter.shape[y]): gridYCenter[:, j, :] = float(j)*cellwidth_y + cellwidth_y/2.0 + \ min_y + grid.lower_bounds[ESMF.StaggerLoc.CENTER_VCENTER][y] * cellwidth_y # last line is the pet specific starting point for this stagger and dim for k in xrange(gridZCenter.shape[z]): gridZCenter[:, :, k] = float(k)*cellwidth_z + cellwidth_z/2.0 + \ min_z + grid.lower_bounds[ESMF.StaggerLoc.CENTER_VCENTER][z] * cellwidth_z # last line is the pet specific starting point for this stagger and dim if domask: mask = grid.add_item(ESMF.GridItem.MASK) mask[:] = 1 mask[np.where((1.75 < gridXCenter.data < 2.25) & (1.75 < gridYCenter.data < 2.25) & (1.75 < gridZCenter.data < 2.25))] = 0 if doarea: grid.add_item(ESMF.GridItem.AREA) area = grid.get_item(ESMF.GridItem.AREA) area[:] = 5.0 return grid
def grid_create_periodic(bounds, domask=False): ''' PRECONDITIONS: 'bounds' contains the number of indices required for the first two dimensions of a periodic Grid. 'domask' is a boolean value that gives the option to put a mask on this Grid.\n POSTCONDITIONS: A periodic Grid has been created.\n RETURN VALUES: \n Grid :: grid \n ''' nx = float(bounds[0]) ny = float(bounds[1]) dx = 360.0/nx dy = 180.0/ny DEG2RAD = 3.141592653589793/180.0 max_index = np.array([nx,ny]) staggerLocs = [ESMF.StaggerLoc.CORNER, ESMF.StaggerLoc.CENTER] grid = ESMF.Grid(max_index, num_peri_dims=1, staggerloc=staggerLocs) # VM vm = ESMF.ESMP_VMGetGlobal() localPet, petCount = ESMF.ESMP_VMGet(vm) # get the coordinate pointers and set the coordinates [x,y] = [0, 1] gridXCorner = grid.get_coords(x, staggerloc=ESMF.StaggerLoc.CORNER) gridYCorner = grid.get_coords(y, staggerloc=ESMF.StaggerLoc.CORNER) # make an array that holds indices from lower_bounds to upper_bounds bnd2indX = np.arange(grid.lower_bounds[ESMF.StaggerLoc.CORNER][x], grid.upper_bounds[ESMF.StaggerLoc.CORNER][x], 1) bnd2indY = np.arange(grid.lower_bounds[ESMF.StaggerLoc.CORNER][y], grid.upper_bounds[ESMF.StaggerLoc.CORNER][y], 1) for i in xrange(gridXCorner.shape[x]): gridXCorner[i, :] = float(bnd2indX[i])*dx - 180.0 for j in xrange(gridYCorner.shape[y]): gridYCorner[:, j] = float(bnd2indY[j])*dy - 90.0 ## CENTERS # get the coordinate pointers and set the coordinates [x,y] = [0, 1] gridXCenter = grid.get_coords(x, staggerloc=ESMF.StaggerLoc.CENTER) gridYCenter = grid.get_coords(y, staggerloc=ESMF.StaggerLoc.CENTER) # make an array that holds indices from lower_bounds to upper_bounds bnd2indX = np.arange(grid.lower_bounds[ESMF.StaggerLoc.CENTER][x], grid.upper_bounds[ESMF.StaggerLoc.CENTER][x], 1) bnd2indY = np.arange(grid.lower_bounds[ESMF.StaggerLoc.CENTER][y], grid.upper_bounds[ESMF.StaggerLoc.CENTER][y], 1) for i in xrange(gridXCenter.shape[x]): gridXCenter[i, :] = float(bnd2indX[i])*dx + 0.5*dx - 180.0 for j in xrange(gridYCenter.shape[y]): y = (float(bnd2indY[j])*dy - 90.0) yp1 = (float(bnd2indY[j]+1)*dy - 90.0) gridYCenter[:, j] = (y+yp1)/2.0 [x,y] = [0, 1] mask = 0 if domask: # set up the grid mask mask = grid.add_item(ESMF.GridItem.MASK) mask[:] = 1 mask[np.where((175. < gridXCenter.data < 185.) & (-5. < gridYCenter.data < 5.))] = 0 return grid
def mesh_create_5(): ''' PRECONDITIONS: None POSTCONDITIONS: A 5 element Mesh has been created. RETURN VALUES: \n Mesh :: mesh \n 4.0 31 ------ 32 ------ 33 | | 22 / | | 21 | / | | | / 23 | 2.0 21 ------ 22 ------ 23 | | | | 11 | 12 | | | | 0.0 11 ------ 12 ------ 13 0.0 2.0 4.0 Node Ids at corners Element Ids in centers Note: This mesh is not parallel, it can only be used in serial ''' # Two parametric dimensions, and two spatial dimensions mesh = ESMF.Mesh(parametric_dim=2, spatial_dim=2) num_node = 9 num_elem = 5 nodeId = np.array([11,12,13,21,22,23,31,32,33]) nodeCoord = np.array([0.0,0.0, # node 11 2.0,0.0, # node 12 4.0,0.0, # node 13 0.0,2.0, # node 21 2.0,2.0, # node 22 4.0,2.0, # node 23 0.0,4.0, # node 31 2.0,4.0, # node 32 4.0,4.0]) # node 33 nodeOwner = np.zeros(num_node) elemId = np.array([11,12,21,22,23]) elemType=np.array([ESMF.MeshElemType.QUAD, ESMF.MeshElemType.QUAD, ESMF.MeshElemType.QUAD, ESMF.MeshElemType.TRI, ESMF.MeshElemType.TRI]) elemConn=np.array([0,1,4,3, # element 11 1,2,5,4, # element 12 3,4,7,6, # element 21 4,8,7, # element 22 4,5,8]) # element 23 mesh.add_nodes(num_node,nodeId,nodeCoord,nodeOwner) mesh.add_elements(num_elem,elemId,elemType,elemConn) return mesh, nodeCoord, nodeOwner, elemType, elemConn
def create_field(grid_or_mesh, name): ''' PRECONDITIONS: An Grid or Mesh has been created, and 'name' is a string that will be used to initialize the name of a new Field.\n POSTCONDITIONS: A Field has been created.\n RETURN VALUES: \n Field :: field \n ''' field = ESMF.Field(grid_or_mesh, name) return field
def initialize_field_grid_periodic(field): ''' PRECONDITIONS: A Field has been created as 'field' with a 'grid' where coordinates have been set on both the center and corner stagger locations. \n POSTCONDITIONS: The 'field' has been initialized to an analytic field.\n RETURN VALUES: \n Field :: field \n ''' DEG2RAD = 3.141592653589793/180.0 # get the coordinate pointers and set the coordinates [x,y] = [0, 1] gridXCoord = field.grid.get_coords(x, ESMF.StaggerLoc.CENTER) gridYCoord = field.grid.get_coords(y, ESMF.StaggerLoc.CENTER) field.data[:] = 2.0 + np.cos(DEG2RAD*gridXCoord)**2 * \ np.cos(2.0*DEG2RAD*(90.0 - gridYCoord)) return field
def run_regridding(srcfield, dstfield, srcfracfield, dstfracfield): ''' PRECONDITIONS: Two Fields have been created and a regridding operation is desired from 'srcfield' to 'dstfield'. The 'srcfracfield' and 'dstfractfield' are Fields created to hold the fractions of the source and destination fields which contribute to conservative regridding.\n POSTCONDITIONS: A regridding operation has set the data on 'dstfield', 'srcfracfield', and 'dstfracfield'.\n RETURN VALUES: \n Field :: dstfield \n Field :: srcfracfield \n Field :: dstfracfield \n ''' # call the regridding functions regridSrc2Dst = ESMF.Regrid(srcfield, dstfield, \ regrid_method=ESMF.RegridMethod.CONSERVE, \ unmapped_action=ESMF.UnmappedAction.ERROR, \ src_frac_field=srcfracfield, \ dst_frac_field=dstfracfield) dstfield = regridSrc2Dst(srcfield, dstfield) return dstfield, srcfracfield, dstfracfield
def compute_mass_grid(valuefield, areafield, dofrac=False, fracfield=None, uninitval=422397696.): ''' PRECONDITIONS: Two Fields have been created and initialized. 'valuefield' contains data values of a field built on the cells of a grid, 'areafield' contains the areas associated with the grid cells, and 'fracfield' contains the fractions of each cell which contributed to a regridding operation involving 'valuefield. 'dofrac' is a boolean value that gives the option to not use the 'fracfield'.\n POSTCONDITIONS: The mass of the data field is computed.\n RETURN VALUES: float :: mass \n ''' mass = 0.0 areafield.get_area() ind = np.where(valuefield.data != uninitval) if dofrac: mass = np.sum(areafield.data[ind] * valuefield.data[ind] * fracfield.data[ind]) else: mass = np.sum(areafield.data[ind] * valuefield.data[ind]) return mass
[1] Khoei S.A. Gharehbaghi A, R. The superconvergent patch recovery technique and data transfer operators in 3d plasticity problems. Finite Elements in Analysis and Design, 43(8), 2007.
[2] K.C. Hung H. Gu, Z. Zong. A modified superconvergent patch recovery method and its application to large deformation problems. Finite Elements in Analysis and Design, 40(5-6), 2004.
[3] D. Ramshaw. Conservative rezoning algorithm for generalized two-dimensional meshes. Journal of Computational Physics, 59, 1985.
[4] Jones, P.W. SCRIP: A Spherical Coordinate Remapping and Interpolation Package. http://www.acl.lanl.gov/climate/software/SCRIP/. Los Alamos National Laboratory Software Release LACC 98-45.