Quad/octree-level classes

The quadtree and octree level classes are used to create, refine and manipulate the semi-structured mesh in parallel. The quadrants or octants are stored in a distributed manner across processors. The following section describes the functions for the forest of octree operations. Analogous functionality is defined for the forest of quadtree object as well. The first step in creating a OctForest object is typically creating a Topology which defines a model with mapped 2D and 3D geometry that consists entirely of non-degenerate quadrilateral surfaces with four non-degenerate edges and hexahedral volumes with 6 enclosing surfaces. This type of object can be created using the constructor as follows:

  • class tmr.TMR.Topology

    The main topology class that contains the objects used to build the underlying mesh.

    This class takes in a general Model, but there are additional requirements that are placed on the model to create a proper topology object. These requirements are as follows:

    1. No edge can degenerate to a vertex.

    2. No face can degenerate to an edge or vertex.

    3. All faces must be surrounded by a single edge loop with 4

    4. All volumes must contain 6 non-degenerate faces that are ordered in coordinate ordering as shown below. Furthermore, all volumes must be of type TFIVolume.

The Topology object defines additional functionality that is generally required only for lower-level operations. Once the topology object has been created, the OctForest can be initialized. It has the following functions:

  • class tmr.TMR.OctForest

    This class defines a forest of octrees. The octrees within the forest can be distributed across processors. The connectivity between octrees is defined on all processors by setting a octree to node connectivity.

    The octrees can be redistributed across processors by using the repartition function. This destroys the nodes that may have been created (but can easily be recomputed).

    balance(self, btype)

    Balance all the elements in the mesh to achieve a 2-to-1 balance

    Parameters

    btype (int) – Indicates whether or not to balance across octant corners

    coarsen(self)

    Create a new forest object by coarsening all the elements within the mesh by one level, if possible. Does not create new nodes.

    Returns

    The coarsened OctForest

    Return type

    OctForest

    createInterpolation(self, forest, vec)

    Create an interpolation object between two octrees that share a common topology.

    Parameters
    • forest (OctForest) – The octree forest that has the common topology

    • vec (VecInterp) – The interpolation operator

    createNodes(self)

    Create and order all the nodes within the mesh

    createRandomTrees(self, nrand=10, min_lev=0, max_lev=8)

    Create a set of trees with random levels of refinement. This is useful for testing purposes.

    Parameters
    • nrand (int) – Number of random octants to create

    • min_lev (int) – Minimum octant refinement level

    • max_lev (int) – Maximum octant refinement level

    createTrees(self, depth=0)

    Create all the octrees in the mesh with the given level of refinement. Be careful, the number of elements created scales with 8**depth!

    Parameters

    depth (int) – Octree of refinement level depth for all trees.

    duplicate(self)

    Duplicate a forest by copying connectivity, topology and elements (but not duplicating nodes)

    Returns

    A duplicate of the OctForest object

    Return type

    OctForest

    getDepNodeConn(self)

    Return the dependent node connectivity, weights and number of dependent nodes from the octree object

    Returns

    Array of dependent nodes, Array of connectivity of dependent nodes Array of weights associated with the dependent nodes

    Return type

    ptr (np.ndarray), conn (np.ndarray), weight (np.ndarray)

    getInterpType(self)

    Get the element interpolation type

    Returns

    The element interpolation type

    Return type

    TMRInterpolationType

    getMeshConn(self)

    Get the portion of the connectivity stored on this processor.

    Returns

    The local part of the connectivity using global node numbers

    Return type

    np.ndarray

    getMeshOrder(self)

    Get the mesh order from the OctForest

    Returns

    Order of the mesh

    Return type

    int

    getNodeRange(self)

    Get the range of nodes each processor owns

    Returns

    An array of ownership ranges

    Return type

    np.ndarray

    getNodesWithName(self, aname)

    Get an array of the node numbers which touch a geometric object with the specified name.

    Parameters

    aname (str) – The name to search for

    Returns

    An array of the local node numbers with name

    Return type

    np.ndarray

    getOctants(self)

    Get an array of the locally owned octants

    Returns

    An array of octants

    Return type

    OctantArray

    getOctsWithName(self, aname)

    Get an array of octants who touch a geometric object with the specified name.

    Parameters

    aname (str) – The name to search for

    Returns

    An array of octants with the specified name

    Return type

    OctantArray

    getPoints(self)

    Get the node locations for all locally owned nodes

    Returns

    An array of node locations

    Return type

    np.ndarray

    getTopology(self)

    Get the Topology object (if any) for this OctForest

    Returns

    Topology object set, None if not set

    Return type

    Topology

    refine(self, refine=None, min_level=0, max_level=MAX_LEVEL)

    Refine the elements in the mesh by the specified number of levels. If a negative number is supplied, coarsen the element.

    Parameters
    • refine (np.ndarray) – Array of integers indicating element refinement

    • min_lev (int) – Minimum octant refinement level

    • max_lev (int) – Maximum octant refinement level

    repartition(self, max_rank=-1)

    Repartition the mesh across processors. This redistributes the elements so that there are an equal, or nearly equal, number of elements on each processor.

    Parameters
    • max_rank (int) – Number of processors to distribute the mesh across.

    • negative, the mesh is distributed across all processors (If) –

    setConnectivity(self, conn)

    Directly set the connectivity in the OctForest. This is in place of the Topology definition.

    Parameters

    conn (np.ndarray) – Connectivity array of shape (n, 8)

    setMeshOrder(self, order, interp=GAUSS_LOBATTO_POINTS)

    Set the order and element nodal pattern of the mesh.

    Parameters
    • order (int) – The number of nodes along an edge

    • interp (TMRInterpolationType) – Type of interpolation to use

    setTopology(self, topo)

    Set the topology of the coarse hexahedral mesh.

    Parameters

    topo (Topology) – Set topo as the underlying topology object

    writeToVTK(self, fname)

    Write the portion of the local forest to a VTK file.

    Parameters

    fname (str) – The file name (unique on each proc)

Similarly, QuadForest can be initialized and has the similar functions:

  • class tmr.TMR.QuadForest

    This class defines a parallel forest of quadrtrees. The connectivity between quadtrees is defined at a global level. The quadrants can easily be redistributed across processors using the repartition() call.

    balance(self, btype)

    Balance all the elements in the mesh to achieve a 2-to-1 balance

    Parameters

    btype (int) – Indicates whether or not to balance across octant corners

    coarsen(self)

    Create a new forest object by coarsening all the elements within the mesh by one level, if possible. Does not create new nodes.

    Returns

    The coarsened QuadForest

    Return type

    QuadForest

    createInterpolation(self, forest, vec)

    Create an interpolation object between two quadtrees that share a common topology.

    Parameters
    • forest (QuadForest) – The quadtree forest that has the common topology

    • vec (VecInterp) – The interpolation operator

    createNodes(self, btype)

    Create and order the nodes in the mesh.

    createRandomTrees(self, nrand=10, min_lev=0, max_lev=8)

    Create a set of trees with random levels of refinement. This is useful for testing purposes.

    Parameters
    • nrand (int) – Number of random quadrants to create

    • min_lev (int) – Minimum quadrant refinement level

    • max_lev (int) – Maximum quadrant refinement level

    createTrees(self, depth=0)

    Create all the quadtrees in the mesh with the given level of refinement. Be careful, the number of elements created scales with 4**depth!

    Parameters

    depth (int) – Level of refinement for all trees.

    duplicate(self)

    Duplicate a forest by copying connectivity, topology and elements (but not duplicating nodes)

    Returns

    A duplicate of the QuadForest object

    Return type

    QuadForest

    getDepNodeConn(self)

    Return the dependent node connectivity, weights and number of dependent nodes from the quadtree object

    Returns

    Array of dependent nodes, Array of connectivity of dependent nodes Array of weights associated with the dependent nodes

    Return type

    ptr (np.ndarray), conn (np.ndarray), weight (np.ndarray)

    getInterpType(self)

    Get the element interpolation type

    Returns

    The element interpolation type

    Return type

    TMRInterpolationType

    getMeshConn(self)

    Get the portion of the connectivity stored on this processor.

    Returns

    The local part of the connectivity using global node numbers

    Return type

    np.ndarray

    getMeshOrder(self)

    Get the mesh order from the QuadForest

    Returns

    Order of the mesh

    Return type

    int

    getNodeRange(self)

    Get the range of nodes each processor owns

    Returns

    An array of ownership ranges

    Return type

    np.ndarray

    getNodesWithName(self, aname)

    Get an array of the node numbers which touch a geometric object with the specified name.

    Parameters

    aname (str) – The name to search for

    Returns

    An array of the local node numbers with name

    Return type

    np.ndarray

    getPoints(self)

    Get the node locations for all locally owned nodes

    Returns

    An array of node locations

    Return type

    np.ndarray

    getQuadrants(self)

    Get an array of the locally owned quadrants

    Returns

    An array of quadrants

    Return type

    QuadrantArray

    getQuadsWithName(self, aname)

    Get an array of quadrants who touch a geometric object with the specified name.

    Parameters

    aname (str) – The name to search for

    Returns

    An array of quadrants with the specified name

    Return type

    QuadrantArray

    getTopology(self)

    Get the Topology object (if any) for this QuadForest

    Returns

    Topology object set, None if not set

    Return type

    Topology

    refine(self, refine=None, min_level=0, max_level=MAX_LEVEL)

    Refine the elements in the mesh by the specified number of levels. If a negative number is supplied, coarsen the element.

    Parameters
    • refine (np.ndarray) – Array of integers indicating element refinement

    • min_lev (int) – Minimum quadrant refinement level

    • max_lev (int) – Maximum quadrant refinement level

    repartition(self)

    Repartition the mesh across processors. This redistributes the elements so that there are an equal, or nearly equal, number of elements on each processor.

    setMeshOrder(self, order, interp=GAUSS_LOBATTO_POINTS)

    Set the order and element nodal pattern of the mesh.

    Parameters
    • order (int) – The number of nodes along an edge

    • interp (TMRInterpolationType) – Type of interpolation to use

    setTopology(self, topo)

    Set the topology of the coarse hexahedral mesh.

    Parameters

    topo (Topology) – Set topo as the underlying topology object

    writeToVTK(self, fname)

    Write the portion of the local forest to a VTK file.

    Parameters

    fname (str) – The file name (unique on each proc)

Typical Usage

The typical usage for a OctForest would consist of the following:

  1. Create the object and call setTopology() to set the super mesh

  2. Create an initial element mesh by calling createTrees()

  3. Create a refined mesh by duplicating and refining the octree forest by calling duplicate() followed by refine(). Note that the length of the integer array passed to refine must be equal to the number of elements in the original mesh.

  4. Balance the mesh and create nodes by calling balance() then createNodes()

  5. Create a sequence of coarser lower-order meshes by repeated calling duplicate() and coarsen()

  6. Create interpolants between meshes for multigrid methods by calling createInterpolation()

See STEP import and meshing for an example of this usage.