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.distmodule contains some operations applying on distributed treesmaia.algo.partmodule 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.
Generalities
Here are some remarks applying to all the functions:
Unless something else specified, functions operate inplace (input tree is modified) and returns
None.The
commargument always refers to the MPI communicator used to create the input tree.Argument
containers_nameis used by some functions transfering data fields (eginterpolate()). Such function operates on a ‘per-container’ basis, meaning that only the requested containers (FlowSolution_t, DiscreteData_t or ZoneSubRegion_t nodes) will be treated. Depending on the function, supported containers can be eitherfull: data exists for all points or elements of all input zones (typically a FlowSolution);
partial: data exists for a susbet of points or elements on some input zones (typically a ZoneSubRegion).
Expected value for
containers_nameis a list ofstror the shortcut'ALL', in which case the fonction selects all the admissible containers.
Distributed algorithms
The following algorithms apply on maia distributed trees.
Connectivities conversions
- convert_s_to_u(dist_tree, connectivity, comm, subset_loc={})
Performs the destructuration of the input distributed tree.
The element connectivity can be generated:
as polyedric elements, if
connectivityis set toPoly;as standard elements, if
connectivityis set toStandard.
In addition, the output location of the
BC_tandGridConnectivity(1to1)_tsubsets can be specified with thesubset_locdictionnary, using respectively the keysBC_tandGC_t. If nothing is specified, these nodes keep their initialGridLocation_tvalue.Input tree is modified inplace.
Note
Exists also as
convert_s_to_ngon()with connectivity set toPolyand subset location set toFaceCenter(orEdgeCenterfor 2D meshes).- Parameters
dist_tree (CGNSDistTree) – Structured distributed tree
connectivity (str) – Type of elements used to describe the connectivity. Admissible values are
"Poly"and"Standard".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 toNone.
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/'S_twoblocks.yaml', MPI.COMM_WORLD) maia.algo.dist.convert_s_to_u(dist_tree, 'Poly', MPI.COMM_WORLD) for zone in maia.pytree.get_all_Zone_t(dist_tree): 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_sortis True, face based PointList keep their original values.Requirement : the
Element_tnodes appearing in the distributed zones must be ordered according to their dimension (either increasing or decreasing). Tree made of Mixed elements are supported (convert_mixed_to_elements()is called under the hood).- Parameters
dist_tree (CGNSDistTree) – 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_sortis 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)
- convert_ngon_to_elements(dist_tree, comm)
Transform a polyedric (NGon based) connectivity into a standard nodal connectivity.
Tree is modified in place : polyedric element, which are supposed to describe only standard elements (tris, quads, tets, pyras, prisms and hexa) are removed and relevant data (such as PointList) are updated.
- Parameters
dist_tree (CGNSDistTree) – distributed tree with polyedric 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, 'Poly', MPI.COMM_WORLD) maia.algo.dist.convert_ngon_to_elements(dist_tree, MPI.COMM_WORLD) elts = PT.get_nodes_from_label(dist_tree, 'Elements_t') assert [PT.Element.Type(e) for e in elts] == ['QUAD_4', 'HEXA_8']
- 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 (CGNSDistTree) – 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 (CGNSDistTree) – 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)
- reorder_elt_sections_from_dim(dist_tree, reverse=False)
Reorder the Elements_t sections of the input tree according to their dimension.
By default, Elements_t nodes are sorted in increasing dimension order (1D, then 2D, then 3D). Decreasing dimension order (3D, then 2D, then 1D) can be obtained using
reverse=True.Input tree is modified inplace.
- Parameters
dist_tree (CGNSDistTree) – Distributed tree
reverse (bool, optional) – If True, elements of the higher dimension get the lower ElementRange. Defaults to
False.
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) for zone in PT.get_all_Zone_t(dist_tree): assert PT.Zone.elt_ordering_by_dim(zone) != 1 # Elts are not increasing by dim maia.algo.dist.reorder_elt_sections_from_dim(dist_tree) for zone in PT.get_all_Zone_t(dist_tree): assert PT.Zone.elt_ordering_by_dim(zone) == 1 # Now, yes
- concatenate_elt_sections(dist_tree, comm)
Gather the Element_t sections of same ElementType into a single one.
Resulting sections are named after their ElementType. Note that :
Sections of same kind must be contiguous to be gathered. This can be achieved using
reorder_elt_sections_from_dim()function.NGON_n,NFACE_nandMIXEDelement kind are not supported.
Input tree is modified inplace.
- Parameters
dist_tree (CGNSDistTree) – Distributed tree
comm (MPIComm) – MPI communicator
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/'H_elt_and_s.yaml', MPI.COMM_WORLD) is_quad_elt = lambda n : PT.get_label(n) == 'Elements_t' and \ PT.Element.Type(n) == 'QUAD_4' assert len(PT.get_nodes_from_predicate(dist_tree, is_quad_elt)) > 1 # Several QUAD sections maia.algo.dist.concatenate_elt_sections(dist_tree, MPI.COMM_WORLD) assert len(PT.get_nodes_from_predicate(dist_tree, is_quad_elt)) == 1 # Now, only one
- generate_jns_vertex_list(dist_tree, comm, have_isolated_faces=False)
For each 1to1 FaceCenter (EdgeCenter in 2D) matching join found in the distributed tree, create a corresponding 1to1 Vertex matching join.
Input tree is modified inplace: Vertex
GridConnectivity_tnodes are created using#Vtxsuffix (long names are hashed).Only unstructured polyedric meshes are supported.
- Parameters
dist_tree (CGNSDistTree) – Distributed tree
comm (MPIComm) – MPI communicator
have_isolated_faces (bool, optional) – Indicate if original joins includes faces who does not share any vertex 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 = maia.io.file_to_dist_tree(mesh_dir/'S_twoblocks.yaml', MPI.COMM_WORLD) maia.algo.dist.convert_s_to_ngon(dist_tree, 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
- find_ridges(dist_tree, bc_identifiers, comm)
Retrieve the edges delimiting specified BC surfaces of a volumic mesh.
Tree is modified inplace: Elements_t nodes containing resulting edge elements are added in input tree.
Setting groups of surfaces
This function retrieves the edges that delimit groups of BCs, which have to be user-provided through the
bc_identifierslist. Each group can be defined by either:the name of each BCs belonging to the group (list of str);
or a family name gathering the BCs (str).
Note that a given BC surface must not appear in more than one group.
Note
For convenience, the shortcut
bc_identifiers='ALL_BCS'can be used to indicate that each BC constitutes an independant group.- Parameters
dist_tree (CGNSDistTree) – Unstructured distributed tree, starting at Zone_t level or higher.
bc_identifiers (list) – List of BC groups bounded by searched edges (see above)
comm (MPIComm) – MPI communicator
Example
import mpi4py.MPI as 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_ATB_45.yaml', MPI.COMM_WORLD) maia.algo.dist.find_ridges(dist_tree, [['wall'], 'AMONT'], MPI.COMM_WORLD) assert PT.get_node_from_path(dist_tree, 'Base/bump_45/topo_edge') is not None
Geometry transformations
- duplicate_from_periodic_jns(dist_tree, zone_paths, jn_paths_for_dupl, dupl_nb, comm, conformize=False, apply_to_fields=True)
Duplicate a mesh from a transformation defined in its periodic connectivities.
Input tree is modified inplace.
- Parameters
dist_tree (CGNSDistTree) – 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.
dupl_nb (int) – Number of duplications to perform
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 toTrue.
See also
For rotating periodicities, it is also possible to automatically recover a circular (360°) mesh with the function
duplicate_from_rotation_jns_to_360(), which takes the same arguments, excepteddupl_nb.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) left_jns = ['Base/bump_45/ZoneGridConnectivity/matchA'] right_jns = ['Base/bump_45/ZoneGridConnectivity/matchB'] maia.algo.dist.duplicate_from_rotation_jns_to_360(dist_tree, ['Base/bump_45'], (left_jns, right_jns), MPI.COMM_WORLD) assert len(maia.pytree.get_all_Zone_t(dist_tree)) == 45
- duplicate_family_from_periodic_jns(dist_tree, family_name, dupl_nb, comm, **kwargs)
Duplicate zones belonging to the specified family.
This is a shortcut for
duplicate_from_periodic_jns()with autodetection of:Zones to duplicate (every zone belonging to the provided family)
Periodic connectivities to use for duplication. Note that this function will fail if the number of periodic transformation found in the group of zones is not exactly one.
- Parameters
dist_tree (CGNSDistTree) – Input distributed tree
family_name (str) – Name of family gathering the zones to duplicate
dupl_nb (int) – Number of duplications to perform
comm (MPIComm) – MPI communicator
kwargs – See
duplicate_from_periodic_jns()
See also
For rotating periodicities, it is also possible to automatically recover a circular (360°) mesh with the function
duplicate_family_from_rotation_jns_to_360(), which takes the same arguments, excepteddupl_nb.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) maia.algo.dist.duplicate_family_from_periodic_jns(dist_tree, 'ATB', 17, MPI.COMM_WORLD) assert len(maia.pytree.get_all_Zone_t(dist_tree)) == 18
- extrude(dist_tree, extrusion_vector, comm, ksubset_as='GC', dupl_vtx_data=False)
Extrude a 2D mesh in the provided direction.
The resulting mesh will be a 3D mesh with one layer of cells. Existing subsets and containers such as
BC_t,FlowSolution_t,ZoneSubRegion_t, etc. are updated following these rules:EdgeCenter regions (lineic) become FaceCenter regions (surfacic),
CellCenter regions (surfacic) become CellCenter regions (volumic),
Vertex regions remain the same if
dupl_vtx_datais False. Otherwise, they are extended with the corresponding extruded vertices, the fields values beeing simply duplicated. This choice does not apply to BC and GC vertex subsets, which are always extended.
In addition, new FaceCenter subsets are created for the initial surface and its extruded counterpart. These subsets can be created as periodic
GridConnectivity_tnodes (usingksubset_as == 'GC') or asBC_tnodes (usingksubset_as == 'BC').Input tree is modified inplace.
- Parameters
dist_tree (CGNSDistTree) – Input 2D distributed tree
extrusion_vector (array of 3 floats) – extrusion axis, which can be any non zero vector
comm (MPIComm) – MPI communicator
ksubset_as (str) – Set kind of surfacic subset created for the initial and extruded planes. Default to
GC.dupl_vtx_data (bool) – Enable duplication of vertex located fields (see above). Default to
False.
Example
from mpi4py import MPI import maia import maia.pytree as PT dist_tree = maia.factory.generate_dist_block(11, 'TRI_3', MPI.COMM_WORLD) maia.algo.dist.extrude(dist_tree, [0,0,0.5], MPI.COMM_WORLD) assert PT.Zone.CellDimension(PT.get_node_from_label(dist_tree, 'Zone_t')) == 3
- merge_zones(dist_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_pathis 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': Subsets having the same name on the original zones (within a same label) produces a unique subset on the output merged zone.'family': Subsets having the same FamilyName on the original zones (within a same label) produces a unique subset on the output merged zone. Subsets without FamilyName fallback to'name'strategy.
Only unstructured polyedric trees are supported, and interfaces between the zones to merge must have a FaceCenter (EdgeCenter in 2D) location.
- Parameters
dist_tree (CGNSDistTree) – 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(dist_tree, family_name, comm, **kwargs)
Merge the zones belonging to the given family into a single one.
See
merge_zones()for full documentation.- Parameters
dist_tree (CGNSDistTree) – Input distributed tree
family_name (str) – Name of the family (read from
(Additional)FamilyName_tnode) 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 theFamily_tnodes 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(dist_tree, comm, **kwargs)
Detect all the zones connected through 1to1 matching jns and merge them.
See
merge_zones()for full documentation.- Parameters
dist_tree (CGNSDistTree) – 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
- remove_degen_faces_from_family(dist_tree, degen_family, comm)
Remove the specified degenerated faces in the input tree.
Degenerated faces are faces whose vertices have distinct ids, but are in fact geometrically reduced to a line or to a single point. This function removes these faces and update the mesh to renumber the other entities. Input tree is modified inplace.
Important
Faces refered by
degen_familymust be degenerated faces, and will be removed anyway.Only U-NGon meshes are managed in this function.
- Parameters
dist_tree (CGNSDistTree) – Input distributed tree, with U-NGon connectivies
degen_family (str) – Name of the family refering to the degenerated faces
comm (MPIComm) – MPI communicator
Example
from mpi4py import MPI import maia from maia.utils.test_utils import sample_mesh_dir dist_tree = maia.io.file_to_dist_tree(sample_mesh_dir/'degen_faces.yaml', MPI.COMM_WORLD) dist_zone = maia.pytree.get_node_from_label(dist_tree, 'Zone_t') #Only one zone assert maia.pytree.Zone.n_face(dist_zone) == 240 maia.algo.dist.remove_degen_faces_from_family(dist_tree, 'DEGEN_AXIS', MPI.COMM_WORLD) assert maia.pytree.Zone.n_face(dist_zone) == 240 - 16 # 16 faces in degenerated family
- 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 (CGNSDistTree) – Input tree
jn_pathes (list of str) – Pathes of the two matching
GridConnectivity_tnodes. Pathes must start from the root of the tree.comm (MPIComm) – MPI communicator
- adapt_mesh_with_feflo(dist_tree, metric, comm, containers_name=[], periodic=False, feflo_opts='', **options)
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
metricargument, 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)
Periodic mesh adaptation
Periodic mesh adaptation is available by activating the
periodicargument. Information from periodic 1to1 GridConnectivity_t nodes in dist_tree will be used to perform mesh adaptation.- Parameters
dist_tree (CGNSDistTree) – 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
containers_name (list of str or
'ALL') – Name of each Vertex located full container to interpolate on the adapted meshperiodic (boolean) – perform periodic mesh adaptation
feflo_opts (str) – Additional arguments passed to Feflo
**options – Additional options (see below)
- Returns
CGNSTree – Adapted mesh (distributed)
The function allows the additional optional parameters:
tmp_dir(str, default to./TMP_adapt_dir) – Absolute or relative path to the directory where meshb files are writtenconstraints(list of str, default to None) : BC names of entities that must not be adapted
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, containers_name=["FlowSolution"], feflo_opts="-c 100 -cmax 100 -p 4")
Interface tools
- connect_1to1_families(dist_tree, families, comm, periodic=None, **options)
Find the matching interface 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 (EdgeCenter in 2D) GridLocation.
If the interface is periodic, the transformation from the first to the second family entities must be specified using the
periodicargument; a dictionnary with keys'translation','rotation_center'and/or'rotation_angle'(in radians) is expected (seemaia.algo.transform_affine()for full description). Missing keys defaults to 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 (resp. edges) remains in their original node which is suffixed by ‘_unmatched’.
This function allows the additional optional parameters:
location– Controls the output GridLocation of the created interfaces. ‘FaceCenter’ (‘EdgeCenter’ in 2D) or ‘Vertex’ are admitted. Defaults to ‘FaceCenter’ or ‘EdgeCenter’, depending of mesh dimension.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 (CGNSDistTree) – Input distributed tree. Only U 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 chosen distribution policy.
Supported policies are:
uniform: each data array is equally reparted over all the processesgather.N: all the data are moved on process Ngather: shortcut forgather.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 (CGNSDistTree) – Distributed tree
policy (str) – distribution policy (see above)
comm (MPIComm) – MPI communicator
Example
from mpi4py import MPI import maia dist_tree = maia.factory.generate_dist_block(21, 'Poly', MPI.COMM_WORLD) maia.algo.dist.redistribute_tree(dist_tree, 'gather.0', MPI.COMM_WORLD)
- concatenate_subsets_from_families(dist_tree, comm, families='*')
For each family, gather the related BC nodes into a single BC.
If the shorcut
'*'is used forfamiliesargument, all the detected FamilyName values in the tree will be used.BCDataSet are concatenated as well, and common metadata (such as FamilyName, Descriptors, etc.) are preserved. Tree is modified inplace and initial BCs are removed after concatenation.
Note
This function add some nodes in resulting BCs to preserve pre-concatenate tree info. Do not delete them if, for any reason, you want to retrieve initial tree (see
deconcatenate_subsets_from_families())If
dist_treehas ZoneSubRegion nodes with BCRegionName related to a concatenated BC, the BCRegionName descriptor will be replaced by the associated PointList.
Warning
For each family-grouped BCs, BCDataSet nodes must have the same tree structure
- Parameters
dist_tree (CGNSDistTree) – Distributed unstructured tree, starting at Zone_t level or higher.
comm (MPIComm) – MPI communicator
families (list of str or '*', optional) – Family names. Default to
"*".
Example
import mpi4py.MPI as 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/'axisym_mesh.yaml', MPI.COMM_WORLD) maia.algo.dist.concatenate_subsets_from_families(dist_tree, MPI.COMM_WORLD, families=['RIDGE']) is_ridge_bc = PT.pred.label_is('BC_t') & PT.pred.belongs_to_family('RIDGE') assert len(PT.get_nodes_from_predicate(dist_tree, is_ridge_bc)) == 1
- deconcatenate_subsets_from_families(dist_tree, comm, families='*')
For each given family, deconcatenate the related BC gathered with the concatenation service.
If the shorcut
'*'is used forfamiliesargument, all the detected FamilyName values in the tree will be used.BCDataSet are deconcatenated as well, and metadata (such as FamilyName, Descriptors, etc.) is preserved on generated BCs. Tree is modified inplace.
Warning
Each family from
familiesargument must lead to unique BC, including the custom node added byconcatenate_subsets_from_families().- Parameters
dist_tree (CGNSDistTree) – Distributed unstructured tree, starting at Zone_t level or higher.
comm (MPIComm) – MPI communicator
families (list of str or '*', optional) – Family names. Default to
"*".
Example
import mpi4py.MPI as 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/'axisym_mesh.yaml', MPI.COMM_WORLD) maia.algo.dist.concatenate_subsets_from_families(dist_tree, MPI.COMM_WORLD, families=['RIDGE']) maia.algo.dist.deconcatenate_subsets_from_families(dist_tree, MPI.COMM_WORLD, families=['RIDGE']) is_ridge_bc = PT.pred.label_is('BC_t') & PT.pred.belongs_to_family('RIDGE') assert len(PT.get_nodes_from_predicate(dist_tree, is_ridge_bc)) == 9
Partitioned algorithms
The following algorithms apply on maia partitioned trees.
Geometric calculations
- 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. BC are considered to be of kind wall if their BCType (or the one of their related family) is one of
'BCWall','BCWallViscous','BCWallViscousHeatFlux'or'BCWallViscousIsothermal'.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 DiscreteData 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 toTrue. Only available when method=cloud.
- Parameters
part_tree (CGNSPartTree) – Input partitioned 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 DiscreteData_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
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=Truepartitioning option.Input field for isosurface computation must be located at vertices.
This function requires ParaDiGMa access.
Note
If
elt_typeis 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.Partial or full containers can be transfered on the output isosurface tree.
Once created, additional fields can be exchanged from volumic tree to isosurface tree using
_exchange_field(part_tree, iso_part_tree, containers_name, comm).
- Parameters
part_tree (CGNSPartTree) – 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 or
'ALL') – Name of each container node to transfer on the output isosurface tree.**options – Options related to plane extraction.
- Returns
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 areTRI_3, QUAD_4, NGON_n. Defaults toTRI_3.graph_part_tool(str) – Controls the isosurface partitioning tool. Admissible values arehilbert, parmetis, ptscotch.hilbertmay produce unbalanced partitions for some configurations. Defaults toptscotch.
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, preserve_orientation=True) 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 (CGNSPartTree) – Partitioned tree to slice. Only U-NGon connectivities are managed.
plane_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 or
'ALL') – Name of each container node to transfer on the output slice tree.**options – Options related to plane extraction (see
iso_surface()).
- Returns
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, preserve_orientation=True) 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 (CGNSPartTree) – Partitioned tree to slice. Only U-NGon connectivities are managed.
sphere_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 or
'ALL') – Name of each container node to transfer on the output slice tree.**options – Options related to plane extraction (see
iso_surface()).
- Returns
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, transfer_dataset=True, 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.
Data fields existing in the volumic mesh can be transfered to the extracted mesh by two ways:
if
transfer_datasetis set toTrue, fields found under the ZoneSubRegion are transfered on the extracted mesh, where they are stored in a FlowSolution_t container since they cover all cells (or vertices).Other full or partial containers are transfered if their name is requested in the
containers_namelist. Their dimension must be at most equal to the one of the extracted mesh.
- Parameters
part_tree (CGNSPartTree) – Partitioned tree from which extraction is computed. U-Elts connectivities are not managed.
zsr_name (str) – Name of the ZoneSubRegion_t node
comm (MPIComm) – MPI communicator
transfer_dataset (bool) – Transfer (or not) fields stored in ZSR to the extracted mesh (default to
True)containers_name (list of str or
'ALL') – Name of each container node to transfer on the output extracted tree.**options – Options related to the extraction.
- Returns
CGNSTree – Extracted submesh (partitioned)
Extraction can be controled by the optional kwargs (only for U meshes):
equilibrate(bool) – IfFalse, the extracted entities remains on their original rank, which simplifies data exchanges but leads to poor load balancing. Default isTrue.graph_part_tool(str) – Partitioning tool used to balance the extracted zones (ifequilibrate=True) Admissible values arehilbert, parmetis, ptscotch. Note that vertex-located extractions require hilbert partitioning. Default ishilbert.
Important
Input tree must have a U-NGon or Structured connectivity
Partitions must come from a single initial domain on input tree.
See also
create_extractor_from_zsr()takes the same parameters, exceptedcontainers_nameandtransfer_dataset, and returns an Extractor object which can be used to exchange containers more than once through itsExtractor.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_namebecomesbc_nameand optionaltransfer_datasetargument allows to transfer BCDataSet (without PointList or PointRange) from BC to the extracted mesh (default toTrue).See also
create_extractor_from_bc_name()takes the same parameters, exceptedcontainers_nameandtransfer_dataset, and returns an Extractor object which can be used to exchange containers more than once through itsExtractor.exchange_fields(container_name)method.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
- extract_part_from_family(part_tree, family_name, comm, transfer_dataset=True, containers_name=[], **options)
Extract the submesh defined by the provided family name from the input volumic partitioned tree.
Family related nodes can be labelled either as BC_t or ZoneSubRegion_t, but their GridLocation must have the same value. They generate a merged output on the resulting extracted tree.
Behaviour and arguments of this function are similar to those of
extract_part_from_zsr().Warning
Only U-NGon meshes are managed in this function.
See also
create_extractor_from_family()takes the same parameters, exceptedcontainers_nameandtransfer_dataset, and returns an Extractor object which can be used to exchange containers more than once through itsExtractor.exchange_fields(container_name)method.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_family(part_tree, \ 'WALL', MPI.COMM_WORLD, containers_name=["WallDistance"]) assert maia.pytree.get_node_from_name(extracted_bc, "WallDistance") is not None
Interpolations
- centers_to_nodes(part_tree, comm, containers_name=[], **options)
Create Vertex located fields from CellCenter located fields.
This transformation is performed for all the fields found under the requested container(s), which must be CellCenter located full containers. Input tree is modified inplace: Vertex containers are created using
#Vtxsuffix.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
part_tree (CGNSPartTree) – Partionned tree
comm (MPIComm) – MPI communicator
containers_name (list of str or
'ALL') – Name of each container node to transfer.**options – Options related to interpolation, see above.
See also
A
CenterToNodeobject can be instanciated with the same parameters, excludingcontainers_name, and then be used to move containers more than once with itsmove_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 fields located at Cells (output container is Geometry_3d) maia.algo.compute_elements_center(part_tree, 3) maia.algo.part.centers_to_nodes(part_tree, comm, ['Geometry_3d']) for part in PT.get_all_Zone_t(part_tree): vtx_sol = PT.get_node_from_name(part, 'Geometry_3d#Vtx') assert PT.Container.GridLocation(vtx_sol) == 'Vertex'
- nodes_to_centers(part_tree, comm, containers_name=[], **options)
Create CellCenter located fields from Vertex located fields.
This transformation is performed for all the fields found under the requested container(s), which must be vertex located full containers. Input tree is modified inplace: CellCenter containers are created using
#Cellsuffix.Interpolation is based on Inverse Distance Weighting (IDW) method: each vertex contributes to the cell value 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.
- Parameters
part_tree (CGNSPartTree) – Partionned tree
comm (MPIComm) – MPI communicator
containers_name (list of str or
'ALL') – Name of each container node to transfer.**options – Options related to interpolation, see above.
See also
A
NodeToCenterobject can be instanciated with the same parameters, excludingcontainers_name, and then be used to move containers more than once with itsmove_fields(container_name)method.Example
import mpi4py import maia import maia.pytree as PT comm = mpi4py.MPI.COMM_WORLD dist_tree = maia.factory.generate_dist_sphere(3, 'TETRA_4', comm) part_tree = maia.factory.partition_dist_tree(dist_tree, comm) # Init a FlowSolution located at Nodes for part in PT.get_all_Zone_t(part_tree): cx, cy, cz = PT.Zone.coordinates(part) fields = {'cX': cx, 'cY': cy, 'cZ': cz} PT.new_FlowSolution('FSol', loc='Vertex', fields=fields, parent=part) maia.algo.part.nodes_to_centers(part_tree, comm, ['FSol']) for part in PT.get_all_Zone_t(part_tree): cell_sol = PT.get_node_from_name(part, 'FSol#Cell') assert PT.Container.GridLocation(cell_sol) == 'CellCenter'
Generic algorithms
The following algorithms apply on maia distributed or partitioned trees
Geometry transformations
- transform_affine(t, rotation_center=None, rotation_angle=None, translation=None, apply_to_fields=True, positional_fields=['Coordinate'])
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 vectors and \(R\) is the rotation matrix. The rotation matrix is computed from
rotation_angle, whose kind depends on physical dimension of the mesh:if
phy_dim == 3, it must be a vector of 3 floats, storing the Euler rotation angles \(\alpha, \beta \text{ and } \gamma\); \(R\) is then the combination ofintrinsic elemental rotations \(X_\alpha, Y^{\prime}_\beta, Z^{\prime\prime}_\gamma\), or, equivalently,
extrinsic elemental rotations \(Z_\gamma, Y_\beta, X_\alpha\).
if
phy_dim == 2, a scalar float \(\theta\) is expected, defining the rotation angle in the XY plane.
Input tree is modified inplace.
- Parameters
t (CGNSTree) – Tree 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 part of the transformation to all the vectorial fields (DataArray_t) found in the input tree. Defaults toTrue.positional_fields (list of str, optional) – If
apply_to_fieldsisTrue, add the translation part for these specific vectorial fields. Defaults to['Coordinate'].
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])
- scale_mesh(t, s=1.0)
Rescale the GridCoordinates of the input mesh.
Input zone(s) can be either structured or unstructured, but must have cartesian coordinates. Transformation is defined by
\[\tilde v = S \cdot v\]where S is the scaling matrix. Input tree is modified inplace.
- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
s (float or array of float) – Scaling factor in each physical dimension. Scalars automatically extend to uniform array.
Example
from mpi4py import MPI import maia dist_tree = maia.factory.generate_dist_block(10, 'Poly', MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(dist_tree, 'CoordinateX')[1].max() <= 1. maia.algo.scale_mesh(dist_tree, [3.0, 2.0, 1.0]) assert maia.pytree.get_node_from_name(dist_tree, 'CoordinateX')[1].max() <= 3.
- cartesian_to_cylindrical(t, axis, comm=None, apply_to_fields=True)
Convert the input tree into a cylindrical coordinate system.
Input zone(s) in the tree can be either structured or unstructured, but must have cartesian coordinates. The revolution axis to be used for the cylindrical coordinate system must be specified thought the
axisargument, and will always be denotedZin the cylindrical system.Input tree is modified inplace; suffixes
R,ThetaandZare used for coordinates and vectorial fields (note that cyl. Z axis and cart. Z axis may differs).- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
axis (array of 3 floats) – Revolution axis, which can by any non zero vector
comm (MPIComm) – MPI communicator, mandatory only for distributed trees
apply_to_fields (bool) – If True, apply the transformation to the vectorial fields found under the following nodes :
FlowSolution_t,DiscreteData_t,ZoneSubRegion_t,BCDataset_t. Defaults toTrue.
Example
import mpi4py.MPI as 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) maia.algo.cartesian_to_cylindrical(dist_tree, (1,0,0), MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(dist_tree, 'CoordinateR') is not None
- cylindrical_to_cartesian(t, axis, comm=None, apply_to_fields=True)
Convert the input tree into a cartesian coordinate system.
Input zone(s) in the tree can be either structured or unstructured, but must have cylindrical coordinates. The expression of the revolution axis (ie the
Zaxis) of the cylindrical coordinate system in the cartesian basis must be specified thought theaxisargument.Input tree is modified inplace; suffixes
R,ThetaandZmust be used for coordinates and vectorial fields (note that cyl. Z axis and cart. Z axis may differs).- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
axis (array of 3 floats) – Revolution axis, which can by any non zero vector
comm (MPIComm) – MPI communicator, mandatory only for distributed trees
apply_to_fields (bool) – If True, apply the transformation to the vectorial fields found under the following nodes :
FlowSolution_t,DiscreteData_t,ZoneSubRegion_t,BCDataset_t. Defaults toTrue.
Example
import mpi4py.MPI as 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_ATB_45.yaml', MPI.COMM_WORLD) maia.algo.cartesian_to_cylindrical(dist_tree, (1,0,0), MPI.COMM_WORLD) # Create a vector field on cylindrical mesh for zone in PT.get_nodes_from_label(dist_tree, 'Zone_t'): cr, ctheta, cz = PT.Zone.coordinates(zone) fields= {'VelocityR': cr**2, 'VelocityTheta': ctheta, 'VelocityZ': 0*cz} PT.new_FlowSolution("FlowSolution", loc="Vertex", fields=fields, parent=zone) maia.algo.cylindrical_to_cartesian(dist_tree, (1,0,0), MPI.COMM_WORLD) assert PT.get_node_from_name(dist_tree, 'VelocityX') is not None
Geometric calculations
- localize_points(src_tree, tgt_tree, location, comm, **options)
Localize points between two 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_tcontainer called “Localization” on the target zones. Note that if the source tree is structured, the output gnum is still a scalar index and not a (i,j,k) triplet.Localization can be parametred thought the options kwargs:
loc_tolerance(default = 1E-6) – Geometric tolerance for the method.
Inputs trees can be either distributed or partitioned, but both must be of same kind.
- Parameters
src_tree (CGNSDistTree|CGNSPartTree) – Source tree
tgt_tree (CGNSDistTree|CGNSPartTree) – Target tree
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 tree_src = maia.io.file_to_dist_tree(mesh_dir/'U_ATB_45.yaml', comm) 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(tree_tgt): maia.algo.transform_affine(tgt_zone, rotation_angle=[170*3.14/180,0,0], translation=[0,0,3]) maia.algo.localize_points(tree_src, tree_tgt, 'CellCenter', comm) for tgt_zone in maia.pytree.get_all_Zone_t(tree_tgt): loc_container = PT.get_child_from_name(tgt_zone, 'Localization') assert PT.Container.GridLocation(loc_container) == 'CellCenter'
- find_closest_points(src_tree, tgt_tree, location, comm)
Find the closest points between two 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_tcontainer called “ClosestPoint” on the target zones. The ids of source points refers to cells or vertices depending on the chosen location.Inputs trees can be either distributed or partitioned, but both must be of same kind.
- Parameters
src_tree (CGNSDistTree|CGNSPartTree) – Source tree
tgt_tree (CGNSDistTree|CGNSPartTree) – Target tree
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.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.Container.GridLocation(loc_container) == 'Vertex'
- compute_elements_center(t, dim, comm=None)
Compute the centers of the specified mesh entity.
The mesh entity on which centers are computed must be specified using
dimparameter: values of 1, 2, and 3 correspond respectively to edges, faces and cells. For convenience, the keywordCellCentercan be used to indicate, on each zone, the higher available dimension. The following table summarizes the possibilities. Note that some combinations do not make sense (zones in this situation are skipped).dim=1
dim=2
dim=3
dim=CellCenter
3D mesh
Edges
Faces
Cells
Cells
2D mesh
Edges
Faces
Faces
1D mesh
Edges
Edges
Warning
For structured meshes,
dim = 1is not yet implemented.Centers are computed using a basic average over the vertices of the entity. Cartesian and cylindrical coordinates are supported.
Input tree is modified inplace : results are stored in a
DiscreteData_tcontainer namedGeometry_{1|2|3}d. Note that for unstructured zones described by standard elements, centers are computed only for elements explicitly defined in sections.- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
dim (int or 'CellCenter') – Entity on which centers are computed (see above)
comm (MPIComm) – MPI communicator, mandatory only for distributed trees
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) maia.algo.compute_elements_center(dist_tree, 3, MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(dist_tree, 'Geometry_3d') is not None
- compute_elements_measure(t, dim, comm=None)
Compute the length, area or volume of the specified mesh entity.
As for
compute_elements_center(), the mesh entity on which measures are computed must be specified usingdimparameter: values of 1, 2, and 3 correspond respectively to edges, faces and cells. For convenience, the keywordCellCentercan be used to indicate, on each zone, the higher available dimension. Seecompute_elements_center()for the summarizing table. Note that some combinations do not make sense (zones in this situation are skipped).Warning
For structured meshes,
dim = 1is not yet implemented.Only cartesian coordinates are supported.
Input tree is modified inplace : results are stored in a
DiscreteData_tcontainer namedGeometry_{1|2|3}d. Note that for unstructured zones described by standard elements, measures are computed only for elements explicitly defined in sections.- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
dim (int or 'CellCenter') – Entity on which measures are computed (see above)
comm (MPIComm) – MPI communicator, mandatory only for distributed trees
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.compute_elements_measure(dist_tree, 3, MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(dist_tree, 'Geometry_3d') is not None
- compute_elements_normal(t, comm=None, unitary=False)
Compute the normal vector of the relevant mesh entity.
Unlike
compute_elements_center()andcompute_elements_measure()this function does not make sense for every mesh entity: actually, this function compute the normal vector toface elements if the physical dimension of the mesh is 3,
edge elements if the physical dimension of the mesh is 2.
Consequently, the number of components of the computed normals equals the physical dimension of the mesh.
If
unitaryisTrue, the norm of the computed vectors is one; otherwise, it is equal to the measure of the corresponding entity.Warning
Only cartesian coordinates are supported.
Input tree is modified inplace : results are stored in a
DiscreteData_tcontainer namedGeometry_{1|2|3}d. Note that for unstructured zones described by standard elements, normals are computed only for elements explicitly defined in sections.- Parameters
t (CGNSTree) – Tree starting at Zone_t level or higher
comm (MPIComm) – MPI communicator, mandatory only for distributed trees
unitary (bool, optional) – If
True, normalize the result (default toFalse).
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.compute_elements_normal(dist_tree, MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(dist_tree, 'Geometry_2d') is not None
Interpolations
- interpolate(src_tree, tgt_tree, comm, containers_name, location, **options)
Interpolate fields between two trees.
This function can transfer CellCenter or Vertex located full container. Target tree is modified inplace: the requested FlowSolution_t containers are transfered from the source tree.
Interpolation strategy can be controled thought the options kwargs:
strategy(default = ‘Closest’) – control interpolation method‘Closest’ : Target points use the inverse distance weighting on the
n_closest_ptsource point values.‘Location’ : For
CellCenterfields, target points take the value of the cell in which they are located. ForVertexfields, target points use finite element weights of source cell vertices to compute interpolation. In both cases, unlocated points take the valueNaN.‘LocationAndClosest’ : Use ‘Location’ method and then ‘ClosestPoint’ method for the unlocated points.
n_closest_pt(default = 1) – If strategy is ‘Closest’ or ‘LocationAndClosest’, specify the number of closest points used for interpolation.loc_tolerance(default = 1E-6) – Geometric tolerance for Location method.
Inputs trees can be either distributed or partitioned, but both must be of same kind.
See also
create_interpolator()takes the same parameters (exceptedcontainers_name, which must be replaced bysrc_location), and returns an Interpolator object which can be used to exchange containers more than once through itsInterpolator.exchange_fields(container_name)method.- Parameters
src_tree (CGNSTree) – Source tree
tgt_tree (CGNSTree) – Target tree
comm (MPIComm) – MPI communicator
containers_name (list of str or
'ALL') – Name of each container node 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.interpolate(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.Container.GridLocation(tgt_sol) == 'Vertex'
Connectivities conversions
- 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) – Distributed, Partitioned or Full tree 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) – Distributed, Partitioned or Full tree 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
- edge_pe_to_ngon(t, comm, removePE=False)
Create a NGon node from a Edge node with ParentElements.
Input tree is modified inplace.
- Parameters
t (CGNSTree) – Distributed, Partitioned or Full tree 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_sphere(5, 'NGON_n', MPI.COMM_WORLD) maia.pytree.rm_nodes_from_name(tree, 'NGonElements') maia.algo.edge_pe_to_ngon(tree, MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(tree, 'NGonElements') is not None
- ngon_to_edge_pe(t, comm, remove_NGon=False)
Create a ParentElements node in the EdgeElements node from a NGon node.
Note that EdgeElement is supposed to exist and define all (including internal) edges. This function retrieves the link between these edges and the NGon node.
Input tree is modified inplace.
- Parameters
t (CGNSTree) – Distributed, Partitioned or Full tree starting at Zone_t level or higher.
comm (MPIComm) – MPI communicator, mandatory only for distributed zones
removeNFace (bool, optional) – If True, remove the NGon node. Defaults to False.
Example
from mpi4py import MPI import maia tree = maia.factory.generate_dist_sphere(5, 'NGON_n', MPI.COMM_WORLD) maia.pytree.rm_nodes_from_name(tree, 'ParentElements') maia.algo.ngon_to_edge_pe(tree, MPI.COMM_WORLD) assert maia.pytree.get_node_from_name(tree, 'ParentElements') is not None
Sequential algorithms
The following algorithms apply on regular (full) pytrees. Note that
these compatibility functions are also wrapped in the maia_poly_old_to_new
and maia_poly_new_to_old scripts, see Quick start section.
- poly_new_to_old(full_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
full_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
Example
import maia from maia.utils.test_utils import mesh_dir tree = maia.io.read_tree(mesh_dir/'U_ATB_45.yaml') assert maia.pytree.get_node_from_name(tree, 'ElementStartOffset') is not None maia.algo.seq.poly_new_to_old(tree) assert maia.pytree.get_node_from_name(tree, 'ElementStartOffset') is None
- poly_old_to_new(full_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
full_tree (CGNSTree) – Tree described with old CGNS convention.
Example
import maia from maia.utils.test_utils import mesh_dir tree = maia.io.read_tree(mesh_dir/'U_ATB_45.yaml') maia.algo.seq.poly_new_to_old(tree) assert maia.pytree.get_node_from_name(tree, 'ElementStartOffset') is None maia.algo.seq.poly_old_to_new(tree) assert maia.pytree.get_node_from_name(tree, 'ElementStartOffset') is not None
- enforce_ngon_pe_local(full_tree)
Shift the ParentElements values in order to make it start at 1, as requested by legacy tools.
The tree is modified in place.
- Parameters
full_tree (CGNSTree) – Tree starting at Zone_t level or higher.
Example
import maia from maia.utils.test_utils import mesh_dir import maia.pytree as PT tree = maia.io.read_tree(mesh_dir/'U_ATB_45.yaml') zone = PT.get_node_from_label(tree, 'Zone_t') n_cell = PT.Zone.n_cell(zone) assert PT.get_node_from_name(zone, 'ParentElements')[1].max() > n_cell maia.algo.seq.enforce_ngon_pe_local(tree) assert PT.get_node_from_name(zone, 'ParentElements')[1].max() <= n_cell