Mesh Conversion

Mesh Converter

The command-line tool MpasMeshConverter.x and the corresponding wrapper function mpas_tools.mesh.conversion.convert() are used to convert a dataset describing cell locations, vertex locations, and connectivity between cells and vertices into a valid MPAS mesh following the MPAS mesh specification.

Example call to the command-line tool:

$ planar_hex --nx 4 --ny 4 --dc 10e3 -o base_mesh.nc
$ MpasMeshConverter.x base_mesh.nc mesh.nc

This example uses the planar_hex tool to generate a small, doubly periodic MPAS mesh with 10-km cells, then uses the MPAS mesh converter to make sure the resulting mesh adheres to the mesh specification. MpasMeshConverter.x takes two arguments, the input and output mesh files, and will prompt the user for file names if these arguments are not supplied.

The same example in a python script can be accomplished with:

from mpas_tools.planar_hex import make_planar_hex_mesh
from mpas_tools.mesh.conversion import convert
from mpas_tools.io import write_netcdf

ds = make_planar_hex_mesh(nx=4, ny=4, dc=10e3, nonperiodic_x=False,
                          nonperiodic_y=False)
ds = convert(ds)
write_netcdf(ds, 'mesh.nc')

Regardless of which of these methods is used, the input mesh must define the following dimensions, variables and global attributes (not the dimension sizes are merely examples from the mesh generated in the previous examples):

netcdf mesh {
dimensions:
  nCells = 16 ;
  nVertices = 32 ;
  vertexDegree = 3 ;
variables:
  double xCell(nCells) ;
  double yCell(nCells) ;
  double zCell(nCells) ;
  double xVertex(nVertices) ;
  double yVertex(nVertices) ;
  double zVertex(nVertices) ;
  int cellsOnVertex(nVertices, vertexDegree) ;
  double meshDensity(nCells) ;

// global attributes:
      :on_a_sphere = "NO" ;
      :sphere_radius = 0. ;
      :is_periodic = "YES" ;

The variable meshDensity is required for historical reasons and is passed unchanged to the resulting mesh. It can contain all ones (or zeros), the resolution of the mesh in kilometers, or whatever field the user wishes.

The following global attributes are optional but will be passed on to the resulting mesh:

// global attributes:
      :x_period = 40000. ;
      :y_period = 34641.0161513775 ;
      :history = "Tue May 26 20:58:10 2020: /home/xylar/miniconda3/envs/mpas/bin/planar_hex --nx 4 --ny 4 --dc 10e3 -o base_mesh.nc" ;

The file_id global attribute is also optional and is preserved in the resulting mesh as parent_id. A new file_id (a random hash tag) is generated by the mesh converter.

The resulting dataset has all the dimensions and variables required for an MPAS mesh.

The mesh converter also generates a file called graph.info that is used to create graph partitions from tools like Metis. This file is not stored by default when the python cull function is used but can be specified with the graphInfoFileName argument. The python function also takes a logger that can be used to capture the output that would otherwise go to the screen via stdout and stderr.

Cell Culler

The command-line tool MpasCellCuller.x and the corresponding wrapper function mpas_tools.mesh.conversion.cull() are used to cull cells from a mesh based on the cullCell field in the input dataset and/or the provided masks. The contents of the cullCell field is merge with the mask(s) from a masking dataset and the inverse of the mask(s) from an inverse-masking dataset. Then, a preserve-masking dataset is used to determine where cells should not be culled.

Example call to the command-line tool, assuming you start with a spherical mesh called base_mesh.nc (not the doubly periodic planar mesh from the examples above):

$ merge_features -c natural_earth -b region -n "Land Coverage" -o land.geojson
$ MpasMaskCreator.x base_mesh.nc land.nc -f land.geojson
$ MpasCellCuller.x base_mesh.nc culled_mesh.nc -m land.nc

This example uses the merge_features tool from the geometric_features conda package to get a geojson file containing land coverage. Then, it uses the mask creator (see the next section) to create a mask on the MPAS mesh that is one inside this region and zero outside. Finally, it culls the base mesh to only those cells where the mask is zero (i.e. the mask indicates which cells are to be removed).

The same example in a python script can be accomplished with:

import xarray
from geometric_features import GeometricFeatures
from mpas_tools.mesh.conversion import mask, cull

gf = GeometricFeatures()

fcLandCoverage = gf.read(componentName='natural_earth', objectType='region',
                         featureNames=['Land Coverage'])

dsBaseMesh = xarray.open_dataset('base_mesh.nc')
dsLandMask = mask(dsBaseMesh, fcMask=fcLandCoverage)
dsCulledMesh = conversion.cull(dsBaseMesh, dsMask=dsLandMask)
write_netcdf(dsCulledMesh, 'culled_mesh.nc')

Here is the full usage of MpasCellCuller.x:

MpasCellCuller.x [input_name] [output_name] [[-m/-i/-p] masks_name] [-c]

    input_name:
        This argument specifies the input MPAS mesh.
    output_name:
        This argument specifies the output culled MPAS mesh.
        If not specified, it defaults to culled_mesh.nc, but
        it is required if additional arguments are specified.
    -m/-i/-p:
        These arguments control how a set of masks is used when
        culling a mesh.
        The -m argument applies a mask to cull based on (i.e.
        where the mask is 1, the mesh will be culled).
        The -i argument applies the inverse mask to cull based
        on (i.e. where the mask is 0, the mesh will be
        culled).
        The -p argument forces any marked cells to not be
        culled.
        If this argument is specified, the masks_name argument
        is required
    -c:
        Output the mapping from old to new mesh (cellMap) in
            cellMapForward.txt,
        and output the reverse mapping from new to old mesh in
            cellMapBackward.txt.

Mask Creator

The command-line tool MpasMaskCreator.x and the corresponding wrapper function mpas_tools.mesh.conversion.mask() are used to create a set of region masks either from mask features or from seed points to be used to flood fill a contiguous block of cells.

Examples usage of the mask creator can be found above under the Cell Culler.

Here is the full usage of MpasMaskCreator.x:

MpasMaskCreator.x in_file out_file [ [-f/-s] file.geojson ] [--positive_lon]
    in_file: This argument defines the input file that masks will be created for.
    out_file: This argument defines the file that masks will be written to.
    -s file.geojson: This argument pair defines a set of points (from the geojson point definition)
        that will be used as seed points in a flood fill algorithim. This is useful when trying to remove isolated cells from a mesh.
    -f file.geojson: This argument pair defines a set of geojson features (regions, transects, or points)
        that will be converted into masks / lists.
    --positive_lon: It is unlikely that you want this argument.  In rare cases when using a non-standard geojson
        file where the logitude ranges from 0 to 360 degrees (with the prime meridian at 0 degrees), use this flag.
        If this flag is not set, the logitude range is -180-180 with 0 degrees being the prime meridian, which is the
        case for standar geojson files including all features from the geometric_feature repo.
        The fact that longitudes in the input MPAS mesh range from 0 to 360 is not relevant to this flag,
        as latitude and longitude are recomputed internally from Cartesian coordinates.
        Whether this flag is passed in or not, any longitudes written are in the 0-360 range.