roms_hgrid module
Tools for creating and working with Arakawa C-Grids
This module is a modified copy from [HHT+21].
- class gridtools.grids.roms_hgrid.BoundaryInteractor(x, y=None, beta=None, ax=None, proj=None, **gridgen_options)
Bases:
object
Interactive grid creation
bry = BoundaryClick(x=[], y=[], beta=None, ax=gca(), **gridgen_options)
The initial boundary polygon points (x and y) are counterclockwise, starting in the upper left corner of the boundary.
Key commands:
t : toggle visibility of verticies d : delete a vertex i : insert a vertex at a point on the polygon line
p : set vertex as beta=1 (a Positive turn, marked with green triangle) m : set vertex as beta=1 (a Negative turn, marked with red triangle) z : set vertex as beta=0 (no corner, marked with a black dot)
G : generate grid from the current boundary using gridgen T : toggle visability of the current grid
- bry.dump(bry_file)
Write the current boundary informtion (bry.x, bry.y, bry.beta) to a cPickle file bry_file.
- bry.load(bry_file)
Read in boundary informtion (x, y, beta) from the cPickle file bry_file.
- bry.remove_grid()
Remove gridlines from axes.
- bry.x
the X boundary points
- bry.y
the Y boundary points
- bry.verts
the verticies of the grid
- bry.grd
the CGrid object
- get_xdata()
- get_ydata()
- load_bry(bry_file='bry.pickle')
- remove_grid()
Remove a generated grid from the BoundaryClick figure
- save_bry(bry_file='bry.pickle')
- save_grid(grid_file='grid.pickle')
- property verts
- property x
- property y
- class gridtools.grids.roms_hgrid.CGrid(x_vert, y_vert, x_rho=None, y_rho=None, x_u=None, y_u=None, x_v=None, y_v=None, x_psi=None, y_psi=None, dx=None, dy=None, dndx=None, dmde=None, angle_rho=None)
Bases:
object
Curvilinear Arakawa C-Grid
The basis for the CGrid class are two arrays defining the verticies of the grid in Cartesian (for geographic coordinates, see CGrid_geo). An optional mask may be defined on the cell centers. Other Arakawa C-grid properties, such as the locations of the cell centers (rho-points), cell edges (u and v velocity points), cell widths (dx and dy) and other metrics (angle, dmde, and dndx) are all calculated internally from the vertex points.
Input vertex arrays may be either type np.array or np.ma.MaskedArray. If masked arrays are used, the mask will be a combination of the specified mask (if given) and the masked locations.
Examples
>>> x, y = mgrid[0.0:7.0, 0.0:8.0] >>> x = np.ma.masked_where( (x<3) & (y<3), x) >>> y = np.ma.MaskedArray(y, x.mask) >>> grd = pyroms.grid.CGrid(x, y) >>> print(grd.x_rho) [[-- -- -- 0.5 0.5 0.5 0.5] [-- -- -- 1.5 1.5 1.5 1.5] [-- -- -- 2.5 2.5 2.5 2.5] [3.5 3.5 3.5 3.5 3.5 3.5 3.5] [4.5 4.5 4.5 4.5 4.5 4.5 4.5] [5.5 5.5 5.5 5.5 5.5 5.5 5.5]] >>> print(grd.mask) [[ 0. 0. 0. 1. 1. 1. 1.] [ 0. 0. 0. 1. 1. 1. 1.] [ 0. 0. 0. 1. 1. 1. 1.] [ 1. 1. 1. 1. 1. 1. 1.] [ 1. 1. 1. 1. 1. 1. 1.] [ 1. 1. 1. 1. 1. 1. 1.]]
- calculate_orthogonality()
Calculate orthogonality error in radians
- property mask
Return mask_rho
- mask_polygon(polyverts, mask_value=0.0)
Mask Cartesian points contained within the polygon defined by polyverts
A cell is masked if the cell center (x_rho, y_rho) is within the polygon. Other sub-masks (mask_u, mask_v, and mask_psi) are updated automatically.
mask_value [=0.0] may be specified to alter the value of the mask set within the polygon. E.g., mask_value=1 for water points.
- property mask_psi
Return mask_psi
- property mask_u
Return mask_u
- property mask_v
Return mask_v
- property x
Return x_vert
- property y
Return x_vert
- class gridtools.grids.roms_hgrid.CGrid_geo(lon_vert, lat_vert, proj, use_gcdist=True, ellipse='WGS84', lon_rho=None, lat_rho=None, lon_u=None, lat_u=None, lon_v=None, lat_v=None, lon_psi=None, lat_psi=None, dx=None, dy=None, dndx=None, dmde=None, angle_rho=None)
Bases:
gridtools.grids.roms_hgrid.CGrid
Curvilinear Arakawa C-grid defined in geographic coordinates
For a geographic grid, a projection may be specified, or The default projection for will be defined by the matplotlib.toolkits.Basemap projection:
proj = Basemap(projection=’merc’, resolution=None, lat_ts=0.0)
For a geographic grid, the cell widths are determined by the great circle distances. Angles, however, are defined using the projected coordinates, so a projection that conserves angles must be used. This means typically either Mercator (projection=’merc’) or Lambert Conformal Conic (projection=’lcc’).
- property lat
Shorthand for lat_vert
- property lon
Shorthand for lon_vert
- mask_polygon_geo(mask_value=0.0)
- class gridtools.grids.roms_hgrid.Focus
Bases:
object
Return a container for a sequence of Focus objects
foc = Focus()
The sequence is populated by using the ‘add_focus_x’ and ‘add_focus_y’ methods. These methods define a point (‘xo’ or ‘yo’), around witch to focus, a focusing factor of ‘focus’, and x and y extent of focusing given by Rx or Ry. The region of focusing will be approximately Gausian, and the resolution will be increased by approximately the value of factor.
- foc.add_focus_x(xo, factor=2.0, Rx=0.1)
- foc.add_focus_y(yo, factor=2.0, Ry=0.1)
- Calls to the object return transformed coordinates:
xf, yf = foc(x, y)
- where x and y must be within [0, 1], and are typically a uniform,
- normalized grid. The focused grid will be the result of applying each of
- the focus elements in the sequence they are added to the series.
Examples
>>> foc = pyroms.grid.Focus() >>> foc.add_focus_x(0.2, factor=3.0, Rx=0.2) >>> foc.add_focus_y(0.6, factor=5.0, Ry=0.35)
>>> x, y = np.mgrid[0:1:3j,0:1:3j] >>> xf, yf = foc(x, y)
>>> print(xf) [[ 0. 0. 0. ] [ 0.36594617 0.36594617 0.36594617] [ 1. 1. 1. ]] >>> print(yf) [[ 0. 0.62479833 1. ] [ 0. 0.62479833 1. ] [ 0. 0.62479833 1. ]]
- add_focus_x(xo, factor=2.0, Rx=0.1)
docstring for add_point
- add_focus_y(yo, factor=2.0, Ry=0.1)
docstring for add_point
- class gridtools.grids.roms_hgrid.Gridgen(xbry, ybry, beta, shape, ul_idx=0, focus=None, proj=None, nnodes=14, precision=1e-12, nppe=3, newton=True, thin=True, checksimplepoly=True, verbose=False)
Bases:
gridtools.grids.roms_hgrid.CGrid
docstring for Gridgen
- generate_grid()
- gridtools.grids.roms_hgrid.dist_point_to_segment(p, s0, s1)
Get the distance of a point to a segment.
*p*, *s0*, *s1* are *xy* sequences
This algorithm is from http://geomalgorithms.com/a02-_lines.html
- class gridtools.grids.roms_hgrid.edit_mask_mesh(grd, proj=None, **kwargs)
Bases:
object
Interactive mask editor
edit_mask_mesh(grd, proj)
Edit grd mask. Mask/Unsmask cell by a simple click on the cell. Mask modification are store in mask_change.txt for further use.
- Key commands:
e : toggle between Editing/Viewing mode
- class gridtools.grids.roms_hgrid.edit_mask_mesh_ij(grd, coast=None, **kwargs)
Bases:
object
Interactive mask editor
edit_mask_mesh_ij(grd)
Edit grd mask. Mask/Unsmask cell by a simple click on the cell. Mask modification are store in mask_change.txt for further use.
- Key commands:
e : toggle between Editing/Viewing mode
- class gridtools.grids.roms_hgrid.get_position_from_map(grd, proj=None, **kwargs)
Bases:
object
Get cell index position Interactively
get_position_from_map(grd, proj)
Get index i, j as well as lon, lat coordinates for one cell simply by clicking on the cell.
- Key commands:
i : toggle between Interactive/Viewing mode
- gridtools.grids.roms_hgrid.rho_to_vert(xr, yr, pm, pn, ang)
- gridtools.grids.roms_hgrid.rho_to_vert_geo(lonr, latr, lonp, latp)
- gridtools.grids.roms_hgrid.uvp_masks(rmask)
return u-, v-, and psi-masks based on input rho-mask
- Parameters
rmask (ndarray) – mask at CGrid rho-points
- Returns
(umask, vmask, pmask) – masks at u-, v-, and psi-points
- Return type
ndarrays