Algo module

The maia.algo module provides various algorithms to be applied to one of the two kind of trees defined by Maia:

  • maia.algo.dist module contains some operations applying on distributed trees

  • maia.algo.part module contains some operations applying on partitioned trees

In addition, some algorithms can be applied indistinctly to distributed or partitioned trees. These algorithms are accessible through the maia.algo module.

The maia.algo.seq module contains a few sequential utility algorithms.

Distributed algorithms

The following algorithms applies on maia distributed trees.

Connectivities conversions

convert_s_to_u(disttree_s, connectivity, comm, subset_loc={})

Performs the destructuration of the input dist_tree.

This function copies the GridCoordinate_t and (full) FlowSolution_t nodes, generate a NGon based connectivity and create a PointList for the following subset nodes: BC_t, BCDataSet_t and GridConnectivity1to1_t. In addition, a PointListDonor node is generated for GridConnectivity_t nodes.

Metadata nodes (“FamilyName_t”, “ReferenceState_t”, …) at zone and base level are also reported on the unstructured tree.

Note

Exists also as convert_s_to_ngon() with connectivity set to NGON_n and subset_loc set to FaceCenter.

Parameters
  • disttree_s (CGNSTree) – Structured tree

  • connectivity (str) – Type of elements used to describe the connectivity. Admissible values are "NGON_n" and "HEXA" (not yet implemented).

  • comm (MPIComm) – MPI communicator

  • subset_loc (dict, optional) – Expected output GridLocation for the following subset nodes: BC_t, GC_t. For each label, output location can be a single location value, a list of locations or None to preserve the input value. Defaults to None.

Returns

CGNSTree – Unstructured disttree

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree_s = maia.io.file_to_dist_tree(mesh_dir/'S_twoblocks.yaml', MPI.COMM_WORLD)

dist_tree_u = maia.algo.dist.convert_s_to_u(dist_tree_s, 'NGON_n', MPI.COMM_WORLD)
for zone in maia.pytree.get_all_Zone_t(dist_tree_u):
  assert maia.pytree.Zone.Type(zone) == "Unstructured"
convert_elements_to_ngon(dist_tree, comm, stable_sort=False)

Transform an element based connectivity into a polyedric (NGon based) connectivity.

Tree is modified in place : standard element are removed from the zones and the PointList are updated. If stable_sort is True, face based PointList keep their original values.

Requirement : the Element_t nodes appearing in the distributed zones must be ordered according to their dimension (either increasing or decreasing).

Parameters
  • dist_tree (CGNSTree) – Tree with connectivity described by standard elements

  • comm (MPIComm) – MPI communicator

  • stable_sort (bool, optional) – If True, 2D elements described in the elements section keep their original id. Defaults to False.

Note that stable_sort is an experimental feature that brings the additional constraints:

  • 2D meshes are not supported;

  • 2D sections must have lower ElementRange than 3D sections.

Example

from mpi4py import MPI
import maia
from maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'Uelt_M6Wing.yaml', MPI.COMM_WORLD)
maia.algo.dist.convert_elements_to_ngon(dist_tree, MPI.COMM_WORLD)
ngons_to_elements(t, comm)

Transform a polyedric (NGon) based connectivity into a standard nodal connectivity.

Tree is modified in place : Polyedric element are removed from the zones and Pointlist (under the BC_t nodes) are updated.

Requirement : polygonal elements are supposed to describe only standard elements (ie tris, quads, tets, pyras, prisms and hexas)

WARNING: this function has not been parallelized yet

Parameters
  • disttree (CGNSTree) – Tree with connectivity described by NGons

  • comm (MPIComm) – MPI communicator

convert_elements_to_mixed(dist_tree, comm)

Transform an element based connectivity into a mixed connectivity.

Tree is modified in place : standard elements are removed from the zones. Note that the original ordering of elements is preserved.

Parameters
  • dist_tree (CGNSTree) – Tree with connectivity described by standard elements

  • comm (MPIComm) – MPI communicator

Example

from mpi4py import MPI
import maia
from maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'Uelt_M6Wing.yaml', MPI.COMM_WORLD)
maia.algo.dist.convert_elements_to_mixed(dist_tree, MPI.COMM_WORLD)
convert_mixed_to_elements(dist_tree, comm)

Transform a mixed connectivity into an element based connectivity.

Tree is modified in place : mixed elements are removed from the zones and the PointList are updated.

Parameters
  • dist_tree (CGNSTree) – Tree with connectivity described by mixed elements

  • comm (MPIComm) – MPI communicator

Example

from mpi4py import MPI
import maia
from maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'Uelt_M6Wing.yaml', MPI.COMM_WORLD)
maia.algo.dist.convert_elements_to_mixed(dist_tree, MPI.COMM_WORLD)
maia.algo.dist.convert_mixed_to_elements(dist_tree, MPI.COMM_WORLD)
rearrange_element_sections(dist_tree, comm)

Rearanges Elements_t sections such that for each zone, sections are ordered in ascending dimensions order and there is only one section by ElementType. Sections are renamed based on their ElementType.

The tree is modified in place. The Elements_t nodes are guaranteed to be ordered by ascending ElementRange.

Parameters
  • dist_tree (CGNSTree) – Tree with an element-based connectivity

  • comm (MPIComm) – MPI communicator

Example

from mpi4py import MPI
import maia
import maia.pytree as PT

dist_tree = maia.factory.generate_dist_block(11, 'PYRA_5', MPI.COMM_WORLD)
pyras = PT.get_node_from_name(dist_tree, 'PYRA_5.0')
assert PT.Element.Range(pyras)[0] == 1 #Until now 3D elements are first

maia.algo.dist.rearrange_element_sections(dist_tree, MPI.COMM_WORLD)
tris = PT.get_node_from_name(dist_tree, 'TRI_3') #Now 2D elements are first
assert PT.Element.Range(tris)[0] == 1
generate_jns_vertex_list(dist_tree, comm, have_isolated_faces=False)

For each 1to1 FaceCenter matching join found in the distributed tree, create a corresponding 1to1 Vertex matching join.

Input tree is modified inplace: Vertex GridConnectivity_t nodes are stored under distinct containers named from the original ones, suffixed with #Vtx. Similarly, vertex GC nodes uses the original name suffixed with #Vtx.

Only unstructured-NGon based meshes are supported.

Parameters
  • dist_tree (CGNSTree) – Distributed tree

  • comm (MPIComm) – MPI communicator

  • have_isolated_faces (bool, optional) – Indicate if original joins includes faces who does not share any edge with other external (join) faces. If False, disable the special treatement needed by such faces (better performances, but will fail if isolated faces were actually present). Defaults to False.

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree_s = maia.io.file_to_dist_tree(mesh_dir/'S_twoblocks.yaml', MPI.COMM_WORLD)
dist_tree = maia.algo.dist.convert_s_to_ngon(dist_tree_s, MPI.COMM_WORLD)

maia.algo.dist.generate_jns_vertex_list(dist_tree, MPI.COMM_WORLD)
assert len(maia.pytree.get_nodes_from_name(dist_tree, 'match*#Vtx')) == 2

Geometry transformations

duplicate_from_rotation_jns_to_360(dist_tree, zone_paths, jn_paths_for_dupl, comm, conformize=False, apply_to_fields=False)

Reconstitute a circular mesh from an angular section of the geometry.

Input tree is modified inplace.

Parameters
  • dist_tree (CGNSTree) – Input distributed tree

  • zone_paths (list of str) – List of pathes (BaseName/ZoneName) of the connected zones to duplicate

  • jn_paths_for_dupl (pair of list of str) – (listA, listB) where listA (resp. list B) stores all the pathes of the GridConnectivity nodes defining the first (resp. second) side of a periodic match.

  • comm (MPIComm) – MPI communicator

  • conformize (bool, optional) – If true, ensure that the generated interface vertices have exactly same coordinates (see conformize_jn_pair()). Defaults to False.

  • apply_to_fields (bool, optional) – See maia.algo.transform_affine(). Defaults to False.

merge_zones(tree, zone_paths, comm, output_path=None, subset_merge='name', concatenate_jns=True)

Merge the given zones into a single one.

Input tree is modified inplace : original zones will be removed from the tree and replaced by the merged zone. Merged zone is added with name MergedZone under the first involved Base except if output_path is not None : in this case, the provided path defines the base and zone name of the merged block.

Subsets of the merged block can be reduced thanks to subset_merge parameter:

  • None : no reduction occurs : all subset of all original zones remains on merged zone, with a numbering suffix.

  • 'name' : Subset having the same name on the original zones (within a same label) produces and unique subset on the output merged zone.

Only unstructured-NGon trees are supported, and interfaces between the zones to merge must have a FaceCenter location.

Parameters
  • tree (CGNSTree) – Input distributed tree

  • zone_paths (list of str) – List of path (BaseName/ZoneName) of the zones to merge. Wildcard * are allowed in BaseName and/or ZoneName.

  • comm (MPIComm) – MPI communicator

  • output_path (str, optional) – Path of the output merged block. Defaults to None.

  • subset_merge (str, optional) – Merging strategy for the subsets. Defaults to ‘name’.

  • concatenate_jns (bool, optional) – if True, reduce the multiple 1to1 matching joins related to the merged_zone to a single one. Defaults to True.

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_Naca0012_multizone.yaml', MPI.COMM_WORLD)
assert len(maia.pytree.get_all_Zone_t(dist_tree)) == 3

maia.algo.dist.merge_zones(dist_tree, ["BaseA/blk1", "BaseB/blk2"], MPI.COMM_WORLD)
assert len(maia.pytree.get_all_Zone_t(dist_tree)) == 2
merge_zones_from_family(tree, family_name, comm, **kwargs)

Merge the zones belonging to the given family into a single one.

See merge_zones() for full documentation.

Parameters
  • tree (CGNSTree) – Input distributed tree

  • family_name (str) – Name of the family (read from FamilyName_t node) used to select the zones.

  • comm (MPIComm) – MPI communicator

  • kwargs – any argument of merge_zones(), excepted output_path

See also

Function merge_all_zones_from_families(tree, comm, **kwargs) does this operation for all the Family_t nodes of the input tree.

Example

from mpi4py import MPI
import maia
import maia.pytree as PT
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_Naca0012_multizone.yaml', MPI.COMM_WORLD)

# FamilyName are not included in the mesh
for zone in PT.get_all_Zone_t(dist_tree):
  PT.new_child(zone, 'FamilyName', 'FamilyName_t', 'Naca0012')  

maia.algo.dist.merge_zones_from_family(dist_tree, 'Naca0012', MPI.COMM_WORLD)

zones = PT.get_all_Zone_t(dist_tree)
assert len(zones) == 1 and PT.get_name(zones[0]) == 'naca0012'
merge_connected_zones(tree, comm, **kwargs)

Detect all the zones connected through 1to1 matching jns and merge them.

See merge_zones() for full documentation.

Parameters
  • tree (CGNSTree) – Input distributed tree

  • comm (MPIComm) – MPI communicator

  • kwargs – any argument of merge_zones(), excepted output_path

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_Naca0012_multizone.yaml', MPI.COMM_WORLD)

maia.algo.dist.merge_connected_zones(dist_tree, MPI.COMM_WORLD)
assert len(maia.pytree.get_all_Zone_t(dist_tree)) == 1
conformize_jn_pair(dist_tree, jn_paths, comm)

Ensure that the vertices belonging to the two sides of a 1to1 GridConnectivity have the same coordinates.

Matching join with Vertex or FaceCenter location are admitted. Coordinates of vertices are made equal by computing the arithmetic mean of the two values.

Input tree is modified inplace.

Parameters
  • dist_tree (CGNSTree) – Input tree

  • jn_pathes (list of str) – Pathes of the two matching GridConnectivity_t nodes. Pathes must start from the root of the tree.

  • comm (MPIComm) – MPI communicator

adapt_mesh_with_feflo(dist_tree, metric, comm, container_names=[], feflo_opts='')

Run a mesh adaptation step using Feflo.a software.

Important

  • Feflo.a is an Inria software which must be installed by you and exposed in your $PATH.

  • This API is experimental. It may change in the future.

Input tree must be unstructured and have an element connectivity. Boundary conditions other than Vertex located are managed.

Adapted mesh is returned as an independant distributed tree.

Setting the metric

Metric choice is available through the metric argument, which can take the following values:

  • None : isotropic adaptation is performed

  • str : path (starting a Zone_t level) to a scalar or tensorial vertex located field:

    • If the path leads to a scalar field (e.g. FlowSolution/Pressure), a metric is computed by Feflo from this field.

    • If the path leads to a tensorial field (e.g. FlowSolution/HessMach), we collect its 6 component (named after CGNS tensor conventions) and pass it to Feflo as a user-defined metric.

    FlowSolution FlowSolution_t
    ├───GridLocation GridLocation_t "Vertex"
    ├───Pressure DataArray_t R8 (100,)
    ├───HessMachXX DataArray_t R8 (100,)
    ├───HessMachYY DataArray_t R8 (100,)
    ├───...
    └───HessMachYZ DataArray_t R8 (100,)
    
  • list of 6 str : each string must be a path to a vertex located field representing one component of the user-defined metric tensor (expected order is XX, XY, XZ, YY, YZ, ZZ)

Parameters
  • dist_tree (CGNSTree) – Distributed tree to be adapted. Only U-Elements single zone trees are managed.

  • metric (str or list) – Path(s) to metric fields (see above)

  • comm (MPIComm) – MPI communicator

  • container_names (list of str) – Name of some Vertex located FlowSolution to project on the adapted mesh

  • feflo_opts (str) – Additional arguments passed to Feflo

Returns

CGNSTree – Adapted mesh (distributed)

Warning

Although this function interface is parallel, keep in mind that Feflo.a is a sequential tool. Input tree is thus internally gathered to a single process, which can cause memory issues on large cases.

Example

import mpi4py.MPI as MPI
import maia
import maia.pytree as PT

from maia.algo.dist import adapt_mesh_with_feflo

dist_tree = maia.factory.generate_dist_block(5, 'TETRA_4', MPI.COMM_WORLD)
zone = PT.get_node_from_label(dist_tree, 'Zone_t')

# > Create a metric field
cx, cy, cz = PT.Zone.coordinates(zone)
fields= {'metric' : (cx-0.5)**5+(cy-0.5)**5 - 1}
PT.new_FlowSolution("FlowSolution", loc="Vertex", fields=fields, parent=zone)

# > Adapt mesh according to scalar metric
adpt_dist_tree = adapt_mesh_with_feflo(dist_tree,
                                       "FlowSolution/metric",
                                       MPI.COMM_WORLD,
                                       container_names=["FlowSolution"],
                                       feflo_opts="-c 100 -cmax 100 -p 4")

Interface tools

connect_1to1_families(dist_tree, families, comm, periodic=None, **options)

Find the matching faces between cgns nodes belonging to the two provided families.

For each one of the two families, all the BC_t or GridConnectivity_t nodes related to the family through a FamilyName/AdditionalFamilyName node will be included in the pairing process. These subset must have a Vertex or FaceCenter GridLocation.

If the interface is periodic, the transformation from the first to the second family entities must be specified using the periodic argument; a dictionnary with keys 'translation', 'rotation_center' and/or 'rotation_angle' is expected. Each key maps to a 3-sized numpy array, with missing keys defaulting zero vector.

Input tree is modified inplace : relevant GridConnectivity_t with PointList and PointListDonor data are created. If all the original elements are successfully paired, the original nodes are removed. Otherwise, unmatched faces remains in their original node which is suffixed by ‘_unmatched’.

This function allows the additional optional parameters:

  • location (default = ‘FaceCenter’) – Controls the output GridLocation of the created interfaces. ‘FaceCenter’ or ‘Vertex’ are admitted.

  • tol (default = 1e-2) – Geometric tolerance used to pair two points. Note that for each vertex, this tolerance is relative to the minimal distance to its neighbouring vertices.

Parameters
  • dist_tree (CGNSTree) – Input distributed tree. Only U-NGon connectivities are managed.

  • families (tuple of str) – Name of the two families to connect.

  • comm (MPIComm) – MPI communicator

  • periodic (dic, optional) – Transformation from first to second family if the interface is periodic. None otherwise. Defaults to None.

  • **options – Additional options

Example

from mpi4py import MPI
from numpy  import array, pi
import maia
import maia.pytree as PT
from maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)

# Remove data that should be created
PT.rm_nodes_from_name(dist_tree, 'PointListDonor')
PT.rm_nodes_from_name(dist_tree, 'GridConnectivityProperty')

# Create FamilyName on interface nodes
PT.new_node('FamilyName', 'FamilyName_t', 'Side1', 
        parent=PT.get_node_from_name(dist_tree, 'matchA'))
PT.new_node('FamilyName', 'FamilyName_t', 'Side2', 
        parent=PT.get_node_from_name(dist_tree, 'matchB'))

maia.algo.dist.connect_1to1_families(dist_tree, ('Side1', 'Side2'), MPI.COMM_WORLD,
        periodic={'rotation_angle' : array([-2*pi/45.,0.,0.])})

assert len(PT.get_nodes_from_name(dist_tree, 'PointListDonor')) == 2

Data management

redistribute_tree(dist_tree, policy, comm)

Redistribute the data of the input tree according to the choosen distribution policy.

Supported policies are:

  • uniform : each data array is equally reparted over all the processes

  • gather.N : all the data are moved on process N

  • gather : shortcut for gather.0

In both case, tree structure remains unchanged on all the processes (the tree is still a valid distributed tree). Input is modified inplace.

Parameters
  • dist_tree (CGNSTree) – Distributed tree

  • policy (str) – distribution policy (see above)

  • comm (MPIComm) – MPI communicator

Example

from mpi4py import MPI
import maia

dist_tree_ini = maia.factory.generate_dist_block(21, 'Poly', MPI.COMM_WORLD)
dist_tree_gathered = maia.algo.dist.redistribute_tree(dist_tree_ini, \
    'gather.0', MPI.COMM_WORLD)

Partitioned algorithms

The following algorithms applies on maia partitioned trees.

Geometric calculations

compute_cell_center(zone)

Compute the cell centers of a partitioned zone.

Input zone must have cartesian coordinates recorded under a unique GridCoordinates node. Centers are computed using a basic average over the vertices of the cells.

Parameters

zone (CGNSTree) – Partitionned CGNS Zone

Returns

cell_center (array) – Flat (interlaced) numpy array of cell centers

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

for zone in maia.pytree.iter_all_Zone_t(part_tree):
  cell_center = maia.algo.part.compute_cell_center(zone)
compute_face_center(zone)

Compute the face centers of a partitioned zone.

Input zone must have cartesian coordinates recorded under a unique GridCoordinates node.

Centers are computed using a basic average over the vertices of the faces.

Note

If zone is described with standard elements, centers will be computed for elements explicitly defined in cgns tree.

Parameters

zone (CGNSTree) – Partitionned 2D or 3D U CGNS Zone

Returns

face_center (array) – Flat (interlaced) numpy array of face centers

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

for zone in maia.pytree.iter_all_Zone_t(part_tree):
  face_center = maia.algo.part.compute_face_center(zone)
compute_edge_center(zone)

Compute the edge centers of a partitioned zone.

Input zone must have cartesian coordinates recorded under a unique GridCoordinates node, and a unstructured standard elements connectivity.

Note

If zone is described with standard elements, centers will be computed for elements explicitly defined in cgns tree.

Parameters

zone (CGNSTree) – Partitionned 2D or 3D U-elts CGNS Zone

Returns

face_center (array) – Flat (interlaced) numpy array of edge centers

Example

from mpi4py import MPI
import maia
dist_tree = maia.factory.generate_dist_block(10, "QUAD_4",  MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

for zone in maia.pytree.iter_all_Zone_t(part_tree):
  edge_center = maia.algo.part.geometry.compute_edge_center(zone)
compute_wall_distance(part_tree, comm, point_cloud='CellCenter', out_fs_name='WallDistance', **options)

Compute wall distances and add it in tree.

For each volumic point, compute the distance to the nearest face belonging to a BC of kind wall. Computation can be done using “cloud” or “propagation” method.

Note

Propagation method requires ParaDiGMa access and is only available for unstructured cell centered NGon connectivities grids. In addition, partitions must have been created from a single initial domain with this method.

Tree is modified inplace: computed distance are added in a FlowSolution container whose name can be specified with out_fs_name parameter.

The following optional parameters can be used to control the underlying method:

  • method ({‘cloud’, ‘propagation’}): Choice of the geometric method. Defaults to 'cloud'.

  • perio (bool): Take into account periodic connectivities. Defaults to True. Only available when method=cloud.

Parameters
  • part_tree (CGNSTree) – Input partitionned tree

  • comm (MPIComm) – MPI communicator

  • point_cloud (str, optional) – Points to project on the surface. Can either be one of “CellCenter” or “Vertex” (coordinates are retrieved from the mesh) or the name of a FlowSolution node in which coordinates are stored. Defaults to CellCenter.

  • out_fs_name (str, optional) – Name of the output FlowSolution_t node storing wall distance data.

  • **options – Additional options related to geometric method (see above)

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

maia.algo.part.compute_wall_distance(part_tree, MPI.COMM_WORLD)
assert maia.pytree.get_node_from_name(part_tree, "WallDistance") is not None
localize_points(src_tree, tgt_tree, location, comm, **options)

Localize points between two partitioned trees.

For all the points of the target tree matching the given location, search the cell of the source tree in which it is enclosed. The result, i.e. the gnum & domain number of the source cell (or -1 if the point is not localized), are stored in a DiscreteData_t container called “Localization” on the target zones.

Source tree must be unstructured and have a NGon connectivity.

Localization can be parametred thought the options kwargs:

  • loc_tolerance (default = 1E-6) – Geometric tolerance for the method.

Parameters
  • src_tree (CGNSTree) – Source tree, partitionned. Only U-NGon connectivities are managed.

  • tgt_tree (CGNSTree) – Target tree, partitionned. Structured or U-NGon connectivities are managed.

  • location ({'CellCenter', 'Vertex'}) – Target points to localize

  • comm (MPIComm) – MPI communicator

  • **options – Additional options related to location strategy

Example

import mpi4py
import maia
import maia.pytree as PT
from   maia.utils.test_utils import mesh_dir
comm = mpi4py.MPI.COMM_WORLD

dist_tree_src = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', comm)
dist_tree_tgt = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', comm)
for tgt_zone in maia.pytree.get_all_Zone_t(dist_tree_tgt):
  maia.algo.transform_affine(tgt_zone, rotation_angle=[170*3.14/180,0,0], translation=[0,0,3])
part_tree_src = maia.factory.partition_dist_tree(dist_tree_src, comm)
part_tree_tgt = maia.factory.partition_dist_tree(dist_tree_tgt, comm)

maia.algo.part.localize_points(part_tree_src, part_tree_tgt, 'CellCenter', comm)
for tgt_zone in maia.pytree.get_all_Zone_t(part_tree_tgt):
  loc_container = PT.get_child_from_name(tgt_zone, 'Localization')
  assert PT.Subset.GridLocation(loc_container) == 'CellCenter'
find_closest_points(src_tree, tgt_tree, location, comm)

Find the closest points between two partitioned trees.

For all points of the target tree matching the given location, search the closest point of same location in the source tree. The result, i.e. the gnum & domain number of the source point, are stored in a DiscreteData_t container called “ClosestPoint” on the target zones. The ids of source points refers to cells or vertices depending on the chosen location.

Parameters
  • src_tree (CGNSTree) – Source tree, partitionned

  • tgt_tree (CGNSTree) – Target tree, partitionned

  • location ({'CellCenter', 'Vertex'}) – Entity to use to compute closest points

  • comm (MPIComm) – MPI communicator

Example

import mpi4py
import maia
import maia.pytree as PT
from   maia.utils.test_utils import mesh_dir
comm = mpi4py.MPI.COMM_WORLD

dist_tree_src = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', comm)
dist_tree_tgt = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', comm)
for tgt_zone in maia.pytree.get_all_Zone_t(dist_tree_tgt):
  maia.algo.transform_affine(tgt_zone, rotation_angle=[170*3.14/180,0,0], translation=[0,0,3])
part_tree_src = maia.factory.partition_dist_tree(dist_tree_src, comm)
part_tree_tgt = maia.factory.partition_dist_tree(dist_tree_tgt, comm)

maia.algo.part.find_closest_points(part_tree_src, part_tree_tgt, 'Vertex', comm)
for tgt_zone in maia.pytree.get_all_Zone_t(part_tree_tgt):
  loc_container = PT.get_child_from_name(tgt_zone, 'ClosestPoint')
  assert PT.Subset.GridLocation(loc_container) == 'Vertex'

Mesh extractions

iso_surface(part_tree, iso_field, comm, iso_val=0.0, containers_name=[], **options)

Create an isosurface from the provided field and value on the input partitioned tree.

Isosurface is returned as an independant (2d) partitioned CGNSTree.

Important

  • Input tree must be unstructured and have a ngon connectivity.

  • Input tree must have been partitioned with preserve_orientation=True partitioning option.

  • Input field for isosurface computation must be located at vertices.

  • This function requires ParaDiGMa access.

Note

  • Once created, additional fields can be exchanged from volumic tree to isosurface tree using _exchange_field(part_tree, iso_part_tree, containers_name, comm).

  • If elt_type is set to ‘TRI_3’, boundaries from volumic mesh are extracted as edges on the isosurface (GridConnectivity_t nodes become BC_t nodes) and FaceCenter fields are allowed to be exchanged.

Parameters
  • part_tree (CGNSTree) – Partitioned tree on which isosurf is computed. Only U-NGon connectivities are managed.

  • iso_field (str) – Path (starting at Zone_t level) of the field to use to compute isosurface.

  • comm (MPIComm) – MPI communicator

  • iso_val (float, optional) – Value to use to compute isosurface. Defaults to 0.

  • containers_name (list of str) – List of the names of the FlowSolution_t nodes to transfer on the output isosurface tree.

  • **options – Options related to plane extraction.

Returns

isosurf_tree (CGNSTree) – Surfacic tree (partitioned)

Isosurface can be controled thought the optional kwargs:

  • elt_type (str) – Controls the shape of elements used to describe the isosurface. Admissible values are TRI_3, QUAD_4, NGON_n. Defaults to TRI_3.

  • graph_part_tool (str) – Controls the isosurface partitioning tool. Admissible values are hilbert, parmetis, ptscotch. hilbert may produce unbalanced partitions for some configurations. Defaults to ptscotch.

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir
dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)
maia.algo.part.compute_wall_distance(part_tree, MPI.COMM_WORLD, point_cloud='Vertex')

part_tree_iso = maia.algo.part.iso_surface(part_tree, "WallDistance/TurbulentDistance", iso_val=0.25,\
    containers_name=['WallDistance'], comm=MPI.COMM_WORLD)

assert maia.pytree.get_node_from_name(part_tree_iso, "WallDistance") is not None
plane_slice(part_tree, plane_eq, comm, containers_name=[], **options)

Create a slice from the provided plane equation \(ax + by + cz - d = 0\) on the input partitioned tree.

Slice is returned as an independant (2d) partitioned CGNSTree. See iso_surface() for use restrictions and additional advices.

Parameters
  • part_tree (CGNSTree) – Partitioned tree to slice. Only U-NGon connectivities are managed.

  • sphere_eq (list of float) – List of 4 floats \([a,b,c,d]\) defining the plane equation.

  • comm (MPIComm) – MPI communicator

  • containers_name (list of str) – List of the names of the FlowSolution_t nodes to transfer on the output slice tree.

  • **options – Options related to plane extraction (see iso_surface()).

Returns

slice_tree (CGNSTree) – Surfacic tree (partitioned)

Example

from mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_Naca0012_multizone.yaml', MPI.COMM_WORLD)
maia.algo.dist.merge_connected_zones(dist_tree, MPI.COMM_WORLD) # Isosurf requires single block mesh
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

slice_tree = maia.algo.part.plane_slice(part_tree, [0,0,1,0.5], MPI.COMM_WORLD, elt_type='QUAD_4')
spherical_slice(part_tree, sphere_eq, comm, containers_name=[], **options)

Create a spherical slice from the provided equation \((x-x_0)^2 + (y-y_0)^2 + (z-z_0)^2 = R^2\) on the input partitioned tree.

Slice is returned as an independant (2d) partitioned CGNSTree. See iso_surface() for use restrictions and additional advices.

Parameters
  • part_tree (CGNSTree) – Partitioned tree to slice. Only U-NGon connectivities are managed.

  • plane_eq (list of float) – List of 4 floats \([x_0, y_0, z_0, R]\) defining the sphere equation.

  • comm (MPIComm) – MPI communicator

  • containers_name (list of str) – List of the names of the FlowSolution_t nodes to transfer on the output slice tree.

  • **options – Options related to plane extraction (see iso_surface()).

Returns

slice_tree (CGNSTree) – Surfacic tree (partitioned)

Example

from mpi4py import MPI
import numpy
import maia
import maia.pytree as PT
dist_tree = maia.factory.generate_dist_block(11, 'Poly', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD, preserve_orientation=True)

# Add solution
zone      = PT.get_node_from_label(part_tree, "Zone_t")
vol_rank  = MPI.COMM_WORLD.Get_rank() * numpy.ones(PT.Zone.n_cell(zone))
src_sol   = PT.new_FlowSolution('FlowSolution', loc='CellCenter', fields={'i_rank' : vol_rank}, parent=zone)

slice_tree = maia.algo.part.spherical_slice(part_tree, [0.5,0.5,0.5,0.25], MPI.COMM_WORLD, \
    ["FlowSolution"], elt_type="NGON_n")

assert maia.pytree.get_node_from_name(slice_tree, "FlowSolution") is not None
extract_part_from_zsr(part_tree, zsr_name, comm, containers_name=[], **options)

Extract the submesh defined by the provided ZoneSubRegion from the input volumic partitioned tree.

Dimension of the output mesh is set up accordingly to the GridLocation of the ZoneSubRegion. Submesh is returned as an independant partitioned CGNSTree and includes the relevant connectivities.

In addition, containers specified in containers_name list are transfered to the extracted tree. Containers to be transfered can be either of label FlowSolution_t or ZoneSubRegion_t.

Parameters
  • part_tree (CGNSTree) – Partitioned tree from which extraction is computed. Only U-NGon connectivities are managed.

  • zsr_name (str) – Name of the ZoneSubRegion_t node

  • comm (MPIComm) – MPI communicator

  • containers_name (list of str) – List of the names of the fields containers to transfer on the output extracted tree.

  • **options – Options related to the extraction.

Returns

extracted_tree (CGNSTree) – Extracted submesh (partitioned)

Extraction can be controled by the optional kwargs:

  • graph_part_tool (str) – Partitioning tool used to balance the extracted zones. Admissible values are hilbert, parmetis, ptscotch. Note that vertex-located extractions require hilbert partitioning. Defaults to hilbert.

Important

  • Input tree must be unstructured and have a ngon connectivity.

  • Partitions must come from a single initial domain on input tree.

See also

create_extractor_from_zsr() takes the same parameters, excepted containers_name, and returns an Extractor object which can be used to exchange containers more than once through its Extractor.exchange_fields(container_name) method.

Example

from   mpi4py import MPI
import numpy as np
import maia
import maia.pytree as PT
from   maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

maia.algo.part.compute_wall_distance(part_tree, MPI.COMM_WORLD, point_cloud='Vertex')

# Create a ZoneSubRegion on procs for extracting odd cells
for part_zone in PT.get_all_Zone_t(part_tree):
  ncell       = PT.Zone.n_cell(part_zone)
  start_range = PT.Element.Range(PT.Zone.NFaceNode(part_zone))[0]
  point_list  = np.arange(start_range, start_range+ncell, 2, dtype=np.int32).reshape((1,-1), order='F')
  PT.new_ZoneSubRegion(name='ZoneSubRegion', point_list=point_list, loc='CellCenter', parent=part_zone)

extracted_tree = maia.algo.part.extract_part_from_zsr(part_tree, 'ZoneSubRegion', MPI.COMM_WORLD,
                                                      containers_name=["WallDistance"])

assert maia.pytree.get_node_from_name(extracted_tree, "WallDistance") is not None
extract_part_from_bc_name(part_tree, bc_name, comm, transfer_dataset=True, containers_name=[], **options)

Extract the submesh defined by the provided BC name from the input volumic partitioned tree.

Behaviour and arguments of this function are similar to those of extract_part_from_zsr() (zsr_name becomes bc_name). Optional transfer_dataset argument allows to transfer BCDataSet from BC to the extracted mesh (default to True).

Example

from   mpi4py import MPI
import maia
from   maia.utils.test_utils import mesh_dir

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', MPI.COMM_WORLD)
part_tree = maia.factory.partition_dist_tree(dist_tree, MPI.COMM_WORLD)

maia.algo.part.compute_wall_distance(part_tree, MPI.COMM_WORLD, point_cloud='Vertex')

extracted_bc = maia.algo.part.extract_part_from_bc_name(part_tree, \
               'wall', MPI.COMM_WORLD, containers_name=["WallDistance"])

assert maia.pytree.get_node_from_name(extracted_bc, "WallDistance") is not None

Interpolations

interpolate_from_part_trees(src_tree, tgt_tree, comm, containers_name, location, **options)

Interpolate fields between two partitionned trees.

For now, interpolation is limited to lowest order: target points take the value of the closest point (or their englobing cell, depending of choosed options) in the source mesh. Interpolation strategy can be controled thought the options kwargs:

  • strategy (default = ‘Closest’) – control interpolation method

    • ‘Closest’ : Target points take the value of the closest source cell center.

    • ‘Location’ : Target points take the value of the cell in which they are located. Unlocated points have take a NaN value.

    • ‘LocationAndClosest’ : Use ‘Location’ method and then ‘ClosestPoint’ method for the unlocated points.

  • loc_tolerance (default = 1E-6) – Geometric tolerance for Location method.

Important

  • Source fields must be located at CellCenter.

  • Source tree must be unstructured and have a ngon connectivity.

See also

create_interpolator_from_part_trees() takes the same parameters, excepted containers_name, and returns an Interpolator object which can be used to exchange containers more than once through its Interpolator.exchange_fields(container_name) method.

Parameters
  • src_tree (CGNSTree) – Source tree, partitionned. Only U-NGon connectivities are managed.

  • tgt_tree (CGNSTree) – Target tree, partitionned. Structured or U-NGon connectivities are managed.

  • comm (MPIComm) – MPI communicator

  • containers_name (list of str) – List of the names of the source FlowSolution_t nodes to transfer.

  • location ({'CellCenter', 'Vertex'}) – Expected target location of the fields.

  • **options – Options related to interpolation strategy

Example

import mpi4py
import numpy
import maia
import maia.pytree as PT
comm = mpi4py.MPI.COMM_WORLD

dist_tree_src = maia.factory.generate_dist_block(11, 'Poly', comm)
dist_tree_tgt = maia.factory.generate_dist_block(20, 'Poly', comm)
part_tree_src = maia.factory.partition_dist_tree(dist_tree_src, comm)
part_tree_tgt = maia.factory.partition_dist_tree(dist_tree_tgt, comm)
# Create fake solution
zone = maia.pytree.get_node_from_label(part_tree_src, "Zone_t")
src_sol = maia.pytree.new_FlowSolution('FlowSolution', loc='CellCenter', parent=zone)
PT.new_DataArray("Field", numpy.random.rand(PT.Zone.n_cell(zone)), parent=src_sol)

maia.algo.part.interpolate_from_part_trees(part_tree_src, part_tree_tgt, comm,\
    ['FlowSolution'], 'Vertex')
tgt_sol = PT.get_node_from_name(part_tree_tgt, 'FlowSolution')
assert tgt_sol is not None and PT.Subset.GridLocation(tgt_sol) == 'Vertex'
centers_to_nodes(tree, comm, containers_name=[], **options)

Create Vertex located FlowSolution_t from CellCenter located FlowSolution_t.

Interpolation is based on Inverse Distance Weighting (IDW) method: each cell contributes to each of its vertices with a weight computed from the distance between the cell isobarycenter and the vertice. The method can be tuned with the following kwargs:

  • idw_power (float, default = 1) – Power to which the cell-vertex distance is elevated.

  • cross_domain (bool, default = True) – If True, vertices located at domain interfaces also receive data from the opposite domain cells. This parameter does not apply to internal partitioning interfaces, which are always crossed.

Parameters
  • tree (CGNSTree) – Partionned tree. Only unstructured connectivities are managed.

  • comm (MPIComm) – MPI communicator

  • containers_name (list of str) – List of the names of the FlowSolution_t nodes to transfer.

  • **options – Options related to interpolation, see above.

See also

A CenterToNode object can be instanciated with the same parameters, excluding containers_name, and then be used to move containers more than once with its move_fields(container_name) method.

Example

import mpi4py
import maia
import maia.pytree as PT
from   maia.utils.test_utils import mesh_dir
comm = mpi4py.MPI.COMM_WORLD

dist_tree = maia.io.file_to_dist_tree(mesh_dir/'U_Naca0012_multizone.yaml', comm)
part_tree = maia.factory.partition_dist_tree(dist_tree, comm)

# Init a FlowSolution located at Cells
for part in PT.get_all_Zone_t(part_tree):
  cell_center = maia.algo.part.compute_cell_center(part)
  fields = {'ccX': cell_center[0::3], 'ccY': cell_center[1::3], 'ccZ': cell_center[2::3]}
  PT.new_FlowSolution('FSol', loc='CellCenter', fields=fields, parent=part)

maia.algo.part.centers_to_nodes(part_tree, comm, ['FSol'])

for part in PT.get_all_Zone_t(part_tree):
  vtx_sol = PT.get_node_from_name(part, 'FSol#Vtx')
  assert PT.Subset.GridLocation(vtx_sol) == 'Vertex'

Generic algorithms

The following algorithms applies on maia distributed or partitioned trees

transform_affine(t, rotation_center=array([0., 0., 0.]), rotation_angle=array([0., 0., 0.]), translation=array([0., 0., 0.]), apply_to_fields=False)

Apply the affine transformation to the coordinates of the given zone.

Input zone(s) can be either structured or unstructured, but must have cartesian coordinates. Transformation is defined by

\[\tilde v = R \cdot (v - c) + c + t\]

where c, t are the rotation center and translation vector and R is the rotation matrix.

Input tree is modified inplace.

Parameters
  • t (CGNSTree(s)) – Tree (or sequences of) starting at Zone_t level or higher.

  • rotation_center (array) – center coordinates of the rotation

  • rotation_angler (array) – angles of the rotation

  • translation (array) – translation vector components

  • apply_to_fields (bool, optional) – if True, apply the rotation vector to the vectorial fields found under following nodes : FlowSolution_t, DiscreteData_t, ZoneSubRegion_t, BCDataset_t. Defaults to False.

Example

from mpi4py import MPI
import maia
dist_tree = maia.factory.generate_dist_block(10, 'Poly', MPI.COMM_WORLD)
zone = maia.pytree.get_all_Zone_t(dist_tree)[0]

maia.algo.transform_affine(zone, translation=[3,0,0])
pe_to_nface(t, comm=None, removePE=False)

Create a NFace node from a NGon node with ParentElements.

Input tree is modified inplace.

Parameters
  • t (CGNSTree(s)) – Distributed or Partitioned tree (or sequences of) starting at Zone_t level or higher.

  • comm (MPIComm) – MPI communicator, mandatory only for distributed zones

  • remove_PE (bool, optional) – If True, remove the ParentElements node. Defaults to False.

Example

from mpi4py import MPI
import maia
tree = maia.factory.generate_dist_block(6, 'Poly', MPI.COMM_WORLD)

for zone in maia.pytree.get_all_Zone_t(tree):
  maia.algo.pe_to_nface(zone, MPI.COMM_WORLD)
  assert maia.pytree.get_child_from_name(zone, 'NFaceElements') is not None
nface_to_pe(t, comm=None, removeNFace=False)

Create a ParentElements node in the NGon node from a NFace node.

Input tree is modified inplace.

Parameters
  • t (CGNSTree(s)) – Distributed or Partitioned tree (or sequences of) starting at Zone_t level or higher.

  • comm (MPIComm) – MPI communicator, mandatory only for distributed zones

  • removeNFace (bool, optional) – If True, remove the NFace node. Defaults to False.

Example

from mpi4py import MPI
import maia
tree = maia.factory.generate_dist_block(6, 'NFace_n', MPI.COMM_WORLD)

maia.algo.nface_to_pe(tree, MPI.COMM_WORLD)
assert maia.pytree.get_node_from_name(tree, 'ParentElements') is not None

Sequential algorithms

The following algorithms applies on regular pytrees.

poly_new_to_old(tree, full_onera_compatibility=True)

Transform a tree with polyhedral unstructured connectivity with new CGNS 4.x conventions to old CGNS 3.x conventions.

The tree is modified in place.

Parameters
  • tree (CGNSTree) – Tree described with new CGNS convention.

  • full_onera_compatibility (bool) – if True, shift NFace and ParentElements ids to begin at 1, irrespective of the NGon and NFace ElementRanges, and make the NFace connectivity unsigned

poly_old_to_new(tree)

Transform a tree with polyhedral unstructured connectivity with old CGNS 3.x conventions to new CGNS 4.x conventions.

The tree is modified in place.

This function accepts trees with old ONERA conventions where NFace and ParentElements ids begin at 1, irrespective of the NGon and NFace ElementRanges, and where the NFace connectivity is unsigned. The resulting tree has the correct CGNS/SIDS conventions.

Parameters

tree (CGNSTree) – Tree described with old CGNS convention.

enforce_ngon_pe_local(t)

Shift the ParentElements values in order to make it start at 1, as requested by legacy tools.

The tree is modified in place.

Parameters

t (CGNSTree(s)) – Tree (or sequences of) starting at Zone_t level or higher.