pyHyp API

class pyhyp.pyHyp.pyHyp(*args, **kwargs)[source]

Create the pyHyp object.

Parameters:
commMPI_INTRACOMM

Comm to use. This is used when running in parallel. If not provided, MPI.COMM_WORLD is used by default.

optionsdict

A dictionary containing the the options for pyHyp.

debugbool

Flag used to specify if debugging. This only needs to be set to true when using a symbolic debugger.

freezeEdge(blockID, edge, dstar)[source]

Specify an edge that will be frozen.

Parameters:
blockIDinteger

Index of block IN ONE BASED ORDERING.

edgestr

String specified for edge. One of ‘ilow’, ‘ihigh’, ‘jlow’, ‘jhigh’

dstartfloat

How much these nodes will influence points around it.

freezeFaces(blockIDs, dstar)[source]

Specify one or more faces (blocks) that will be frozen

Parameters:
blockIDsinteger or list

Index of block(s) IN ONE BASED ORDERING.

dstartfloat

How much these nodes will influence points around it.

getSurfaceCoordinates()[source]

Return the surface coordinates on this processor

run()[source]

Run given using the options given

setSurfaceCoordinates(coords)[source]

Set the surface coordinates on this processor

surfaceSmooth(nIter, stepSize, surfFile=None)[source]

Run smoothing iterations on the body surface

Parameters:
nIterint

Number of iterations to run

stepSizefloat

Size of step. Must be < 1. Usually less than 0.1 for stability reasons.

writeCGNS(fileName)[source]

After we have generated a grid, write it out in a properly formatted 1-Cell wide CGNS file suitable for running in SUmb.

writeLayer(fileName, layer=1, meshType='plot3d', partitions=True)[source]

Write a single mesh layer out to a file for visualization or for other purposes.

Parameters:
fileNamestr

Filename to use. Should have .fmt extension for plot3d or .dat for tecplot

layerint

Index of layer to print. Values greater than 1 are only valid if the mesh has already been extruded.

meshTypestr

Type of mesh to write. The two valid arguments are ‘plot3d’ and ‘fe’. The plot3d will write the mesh in the original plot3d format while the FE mesh is the unstructured internal representation.

partitionsbool

This flag which is only used for the ‘fe’ mesh type option will write a separate zone for each partition on each processor. This is useful for visualizing the parallel mesh decomposition.

writeOutput(fileName=None, fileType=None)[source]

This selects the output type based on what is specified in the options

writePlot3D(fileName)[source]

After we have generated a grid, write it out to a plot3d file for the user to look at

class pyhyp.pyHyp.pyHypMulti(comm=None, options=None, commonOptions=None, debug=False, skipList=[])[source]

This is class can be used to run multiple pyHyp cases at once.

The inititalization method will setup, run, and write all the results.

Parameters:
optionsobject

ORDERED dictionary or list of dictionaries. This contains options for the extrusion of all several grids. An example of option dictionary is given below:

options = {
    "epsE": 4.0,
    "epsI": 8.0,
    "outputFile": "corner_hyp.cgns",
    "skip": False,
}

We can set a list of dictionaries as input:

options1 = {
    "epsE": 4.0,
    "epsI": 8.0,
    "outputFile": "corner1_hyp.cgns",
    "skip": False,
}
options2 = "cartesian.cgns"
options3 = {
    "epsE": 2.0,
    "epsI": 4.0,
    "outputFile": "corner2_hyp.cgns",
    "skip": False,
}
options = [options1, options2, options3]

Alternatively, we can set an ORDERED dictionary of dictionaries as input:

from collections import OrderedDict

options = OrderedDict()
options["case1"] = {
    "epsE": 4.0,
    "epsI": 8.0,
    "outputFile": "corner1_hyp.cgns",
    "skip": False,
}
options["block"] = "cartesian.cgns"
options["case2"] = {
    "epsE": 2.0,
    "epsI": 4.0,
    "outputFile": "corner2_hyp.cgns",
    "skip": False,
}

Each element of the list/dictionary will be considered as a different case. One of the elements can be a string specifying a CGNS file that should be combined with the other grids in the end. pyHyp will not do anything with this file except combine it with the generated grids in the corresponding order. These options will overwrite the default options (defined in the pyHyp class) and the common options (another argument of this method). If the user gives a list, this will be converted to a dictionary with integers as keys. Remember this when setting the skip list for unnamed cases.

commomOptionsdict

Dictionary with options that should be applied to all cases in the options dictionary. See the ‘defOpts’ dictionary defined in the pyHyp class to see the available options.

skip_listlist

List containing names of cases that should be skipped.

combineCGNS(combinedFile='combined.cgns', additionalGrids=[], skipList=[], eraseFiles=True)[source]

This will gather all newly generated grids and combine them in a single CGNS file. This only works for CGNS output files.

Parameters:
combinedFilestr

The name of the combined output file.

additionalGridslist

The filenames of any grids that were not generated by the current pyHypMulti object, but should still be included in the combined output file.

skipListlist

The keys of any generated grids that should not be included in the combined output file.

eraseFilesbool

If True, we erase the individual files that are combined.

writeLog()[source]

This will print a log with important information regarding all grids

pyhyp.utils.simpleOCart(inputGrid, dh, hExtra, nExtra, sym, mgcycle, outFile, userOptions=None, xBounds=None, useFarfield=True)[source]

Generates a Cartesian mesh around the provided grid, surrounded by an O-mesh.

Parameters:
inputGridcgnsutils Grid object or str

If a cgnsutils Grid object is provided, we use it as is. Alternatively if a string is provided, we treat it as the name of the nearfield CGNS file to mesh around.

dhfloat or list of float

The target edge length of each cell in the Cartesian part of the mesh. A list of (x, y, z) lengths can be provided to make non-cubic cells. The actual edge lengths will depend on mgcycle.

hExtrafloat

The distance from the Cartesian mesh boundary to the farfield.

nExtraint

The number of layers to extrude the hyperbolic O-mesh.

symstr or list of str

Axis or plane of symmetry. One or more of (‘x’, ‘y’, ‘z’, ‘xmin’, ‘xmax’, ‘ymin’, ‘ymax’, ‘zmin’, ‘zmax’).

mgcycleint or list of int

Number of times mesh should be able to be coarsened for multigrid cycles. A list can be provided for nonuniform (x, y, z) coarsening.

outFilestr

Output file name.

userOptionsdict, optional

Custom pyhyp options to be used with this extrusion. If overset BCs are desired on the outer face, do not set it in this dictionary because after extrusion we overwrite all BCs on the combined grid. See the option useFarfield below. The default value (True) results in farfield BCs on the outer face, and setting it to false results in overset for the far face. Other pyhyp extrusion parameters can be set here.

xBoundsarray (2 x 3), optional

Optional bounding box coordinates desired for the center cartesian grid. The default value can be obtained by: xMin, xMax = grid.getBoundingBox(), and then the xBounds array can be set as xBounds = [xMin, xMax]. This option allows users to modify the bounding box coordinates rather than simply defaulting to the bounding box of the nearfield grid.

useFarfieldbool, optional

Optional flag to control the outermost layer’s BC. Default, True, will result in a farfield outer layer, setting this to False does overset BCs on the outermost layer