Mesh Functions

gmsh.model.mesh.optimizeFunction
gmsh.model.mesh.optimize(method, force = false, niter = 1, dimTags = Tuple{Cint,Cint}[])

Optimize the mesh of the current model using method (empty for default tetrahedral mesh optimizer, "Netgen" for Netgen optimizer, "HighOrder" for direct high-order mesh optimizer, "HighOrderElastic" for high-order elastic smoother, "HighOrderFastCurving" for fast curving algorithm, "Laplace2D" for Laplace smoothing, "Relocate2D" and "Relocate3D" for node relocation). If force is set apply the optimization also to discrete entities. If dimTags is given, only apply the optimizer to the given entities.

source
gmsh.model.mesh.getLastEntityErrorFunction
gmsh.model.mesh.getLastEntityError()

Get the last entities (if any) where a meshing error occurred. Currently only populated by the new 3D meshing algorithms.

Return dimTags.

source
gmsh.model.mesh.getLastNodeErrorFunction
gmsh.model.mesh.getLastNodeError()

Get the last nodes (if any) where a meshing error occurred. Currently only populated by the new 3D meshing algorithms.

Return nodeTags.

source
gmsh.model.mesh.clearFunction
gmsh.model.mesh.clear(dimTags = Tuple{Cint,Cint}[])

Clear the mesh, i.e. delete all the nodes and elements, for the entities dimTags. if dimTags is empty, clear the whole mesh. Note that the mesh of an entity can only be cleared if this entity is not on the boundary of another entity with a non-empty mesh.

source
gmsh.model.mesh.getNodesFunction
gmsh.model.mesh.getNodes(dim = -1, tag = -1, includeBoundary = false, returnParametricCoord = true)

Get the nodes classified on the entity of dimension dim and tag tag. If tag < 0, get the nodes for all entities of dimension dim. If dim and tag are negative, get all the nodes in the mesh. nodeTags contains the node tags (their unique, strictly positive identification numbers). coord is a vector of length 3 times the length of nodeTags that contains the x, y, z coordinates of the nodes, concatenated: [n1x, n1y, n1z, n2x, ...]. If dim >= 0 and returnParamtricCoord is set, parametricCoord contains the parametric coordinates ([u1, u2, ...] or [u1, v1, u2, ...]) of the nodes, if available. The length of parametricCoord can be 0 or dim times the length of nodeTags. If includeBoundary is set, also return the nodes classified on the boundary of the entity (which will be reparametrized on the entity if dim >= 0 in order to compute their parametric coordinates).

Return nodeTags, coord, parametricCoord.

source
gmsh.model.mesh.getNodesByElementTypeFunction
gmsh.model.mesh.getNodesByElementType(elementType, tag = -1, returnParametricCoord = true)

Get the nodes classified on the entity of tag tag, for all the elements of type elementType. The other arguments are treated as in getNodes.

Return nodeTags, coord, parametricCoord.

source
gmsh.model.mesh.getNodeFunction
gmsh.model.mesh.getNode(nodeTag)

Get the coordinates and the parametric coordinates (if any) of the node with tag tag. This function relies on an internal cache (a vector in case of dense node numbering, a map otherwise); for large meshes accessing nodes in bulk is often preferable.

Return coord, parametricCoord.

source
gmsh.model.mesh.setNodeFunction
gmsh.model.mesh.setNode(nodeTag, coord, parametricCoord)

Set the coordinates and the parametric coordinates (if any) of the node with tag tag. This function relies on an internal cache (a vector in case of dense node numbering, a map otherwise); for large meshes accessing nodes in bulk is often preferable.

source
gmsh.model.mesh.getNodesForPhysicalGroupFunction
gmsh.model.mesh.getNodesForPhysicalGroup(dim, tag)

Get the nodes from all the elements belonging to the physical group of dimension dim and tag tag. nodeTags contains the node tags; coord is a vector of length 3 times the length of nodeTags that contains the x, y, z coordinates of the nodes, concatenated: [n1x, n1y, n1z, n2x, ...].

Return nodeTags, coord.

source
gmsh.model.mesh.addNodesFunction
gmsh.model.mesh.addNodes(dim, tag, nodeTags, coord, parametricCoord = Cdouble[])

Add nodes classified on the model entity of dimension dim and tag tag. nodeTags contains the node tags (their unique, strictly positive identification numbers). coord is a vector of length 3 times the length of nodeTags that contains the x, y, z coordinates of the nodes, concatenated: [n1x, n1y, n1z, n2x, ...]. The optional parametricCoord vector contains the parametric coordinates of the nodes, if any. The length of parametricCoord can be 0 or dim times the length of nodeTags. If the nodeTags vector is empty, new tags are automatically assigned to the nodes.

source
gmsh.model.mesh.reclassifyNodesFunction
gmsh.model.mesh.reclassifyNodes()

Reclassify all nodes on their associated model entity, based on the elements. Can be used when importing nodes in bulk (e.g. by associating them all to a single volume), to reclassify them correctly on model surfaces, curves, etc. after the elements have been set.

source
gmsh.model.mesh.relocateNodesFunction
gmsh.model.mesh.relocateNodes(dim = -1, tag = -1)

Relocate the nodes classified on the entity of dimension dim and tag tag using their parametric coordinates. If tag < 0, relocate the nodes for all entities of dimension dim. If dim and tag are negative, relocate all the nodes in the mesh.

source
gmsh.model.mesh.getElementsFunction
gmsh.model.mesh.getElements(dim = -1, tag = -1)

Get the elements classified on the entity of dimension dim and tag tag. If tag < 0, get the elements for all entities of dimension dim. If dim and tag are negative, get all the elements in the mesh. elementTypes contains the MSH types of the elements (e.g. 2 for 3-node triangles: see getElementProperties to obtain the properties for a given element type). elementTags is a vector of the same length as elementTypes; each entry is a vector containing the tags (unique, strictly positive identifiers) of the elements of the corresponding type. nodeTags is also a vector of the same length as elementTypes; each entry is a vector of length equal to the number of elements of the given type times the number N of nodes for this type of element, that contains the node tags of all the elements of the given type, concatenated: [e1n1, e1n2, ..., e1nN, e2n1, ...].

Return elementTypes, elementTags, nodeTags.

source
gmsh.model.mesh.getElementFunction
gmsh.model.mesh.getElement(elementTag)

Get the type and node tags of the element with tag tag. This function relies on an internal cache (a vector in case of dense element numbering, a map otherwise); for large meshes accessing elements in bulk is often preferable.

Return elementType, nodeTags.

source
gmsh.model.mesh.getElementByCoordinatesFunction
gmsh.model.mesh.getElementByCoordinates(x, y, z, dim = -1, strict = false)

Search the mesh for an element located at coordinates (x, y, z). This function performs a search in a spatial octree. If an element is found, return its tag, type and node tags, as well as the local coordinates (u, v, w) within the reference element corresponding to search location. If dim is >= 0, only search for elements of the given dimension. If strict is not set, use a tolerance to find elements near the search location.

Return elementTag, elementType, nodeTags, u, v, w.

source
gmsh.model.mesh.getElementsByCoordinatesFunction
gmsh.model.mesh.getElementsByCoordinates(x, y, z, dim = -1, strict = false)

Search the mesh for element(s) located at coordinates (x, y, z). This function performs a search in a spatial octree. Return the tags of all found elements in elementTags. Additional information about the elements can be accessed through getElement and getLocalCoordinatesInElement. If dim is >= 0, only search for elements of the given dimension. If strict is not set, use a tolerance to find elements near the search location.

Return elementTags.

source
gmsh.model.mesh.getLocalCoordinatesInElementFunction
gmsh.model.mesh.getLocalCoordinatesInElement(elementTag, x, y, z)

Return the local coordinates (u, v, w) within the element elementTag corresponding to the model coordinates (x, y, z). This function relies on an internal cache (a vector in case of dense element numbering, a map otherwise); for large meshes accessing elements in bulk is often preferable.

Return u, v, w.

source
gmsh.model.mesh.getElementTypesFunction
gmsh.model.mesh.getElementTypes(dim = -1, tag = -1)

Get the types of elements in the entity of dimension dim and tag tag. If tag < 0, get the types for all entities of dimension dim. If dim and tag are negative, get all the types in the mesh.

Return elementTypes.

source
gmsh.model.mesh.getElementTypeFunction
gmsh.model.mesh.getElementType(familyName, order, serendip = false)

Return an element type given its family name familyName ("Point", "Line", "Triangle", "Quadrangle", "Tetrahedron", "Pyramid", "Prism", "Hexahedron") and polynomial order order. If serendip is true, return the corresponding serendip element type (element without interior nodes).

Return an integer value.

source
gmsh.model.mesh.getElementPropertiesFunction
gmsh.model.mesh.getElementProperties(elementType)

Get the properties of an element of type elementType: its name (elementName), dimension (dim), order (order), number of nodes (numNodes), local coordinates of the nodes in the reference element (localNodeCoord vector, of length dim times numNodes) and number of primary (first order) nodes (numPrimaryNodes).

Return elementName, dim, order, numNodes, localNodeCoord, numPrimaryNodes.

source
gmsh.model.mesh.getElementsByTypeFunction
gmsh.model.mesh.getElementsByType(elementType, tag = -1, task = 0, numTasks = 1)

Get the elements of type elementType classified on the entity of tag tag. If tag < 0, get the elements for all entities. elementTags is a vector containing the tags (unique, strictly positive identifiers) of the elements of the corresponding type. nodeTags is a vector of length equal to the number of elements of the given type times the number N of nodes for this type of element, that contains the node tags of all the elements of the given type, concatenated: [e1n1, e1n2, ..., e1nN, e2n1, ...]. If numTasks > 1, only compute and return the part of the data indexed by task.

Return elementTags, nodeTags.

source
gmsh.model.mesh.addElementsFunction
gmsh.model.mesh.addElements(dim, tag, elementTypes, elementTags, nodeTags)

Add elements classified on the entity of dimension dim and tag tag. types contains the MSH types of the elements (e.g. 2 for 3-node triangles: see the Gmsh reference manual). elementTags is a vector of the same length as types; each entry is a vector containing the tags (unique, strictly positive identifiers) of the elements of the corresponding type. nodeTags is also a vector of the same length as types; each entry is a vector of length equal to the number of elements of the given type times the number N of nodes per element, that contains the node tags of all the elements of the given type, concatenated: [e1n1, e1n2, ..., e1nN, e2n1, ...].

source
gmsh.model.mesh.addElementsByTypeFunction
gmsh.model.mesh.addElementsByType(tag, elementType, elementTags, nodeTags)

Add elements of type elementType classified on the entity of tag tag. elementTags contains the tags (unique, strictly positive identifiers) of the elements of the corresponding type. nodeTags is a vector of length equal to the number of elements times the number N of nodes per element, that contains the node tags of all the elements, concatenated: [e1n1, e1n2, ..., e1nN, e2n1, ...]. If the elementTag vector is empty, new tags are automatically assigned to the elements.

source
gmsh.model.mesh.getIntegrationPointsFunction
gmsh.model.mesh.getIntegrationPoints(elementType, integrationType)

Get the numerical quadrature information for the given element type elementType and integration rule integrationType (e.g. "Gauss4" for a Gauss quadrature suited for integrating 4th order polynomials). localCoord contains the u, v, w coordinates of the G integration points in the reference element: [g1u, g1v, g1w, ..., gGu, gGv, gGw]. weights contains the associated weights: [g1q, ..., gGq].

Return localCoord, weights.

source
gmsh.model.mesh.getJacobiansFunction
gmsh.model.mesh.getJacobians(elementType, localCoord, tag = -1, task = 0, numTasks = 1)

Get the Jacobians of all the elements of type elementType classified on the entity of tag tag, at the G evaluation points localCoord given as concatenated triplets of coordinates in the reference element [g1u, g1v, g1w, ..., gGu, gGv, gGw]. Data is returned by element, with elements in the same order as in getElements and getElementsByType. jacobians contains for each element the 9 entries of the 3x3 Jacobian matrix at each evaluation point. The matrix is returned by column: [e1g1Jxu, e1g1Jyu, e1g1Jzu, e1g1Jxv, ..., e1g1Jzw, e1g2Jxu, ..., e1gGJzw, e2g1Jxu, ...], with Jxu=dx/du, Jyu=dy/du, etc. determinants contains for each element the determinant of the Jacobian matrix at each evaluation point: [e1g1, e1g2, ... e1gG, e2g1, ...]. coord contains for each element the x, y, z coordinates of the evaluation points. If tag < 0, get the Jacobian data for all entities. If numTasks > 1, only compute and return the part of the data indexed by task.

Return jacobians, determinants, coord.

source
gmsh.model.mesh.getJacobianFunction
gmsh.model.mesh.getJacobian(elementTag, localCoord)

Get the Jacobian for a single element elementTag, at the G evaluation points localCoord given as concatenated triplets of coordinates in the reference element [g1u, g1v, g1w, ..., gGu, gGv, gGw]. jacobians contains the 9 entries of the 3x3 Jacobian matrix at each evaluation point. The matrix is returned by column: [e1g1Jxu, e1g1Jyu, e1g1Jzu, e1g1Jxv, ..., e1g1Jzw, e1g2Jxu, ..., e1gGJzw, e2g1Jxu, ...], with Jxu=dx/du, Jyu=dy/du, etc. determinants contains the determinant of the Jacobian matrix at each evaluation point. coord contains the x, y, z coordinates of the evaluation points. This function relies on an internal cache (a vector in case of dense element numbering, a map otherwise); for large meshes accessing Jacobians in bulk is often preferable.

Return jacobians, determinants, coord.

source
gmsh.model.mesh.getBasisFunctionsFunction
gmsh.model.mesh.getBasisFunctions(elementType, localCoord, functionSpaceType, wantedOrientations = Cint[])

Get the basis functions of the element of type elementType at the evaluation points localCoord (given as concatenated triplets of coordinates in the reference element [g1u, g1v, g1w, ..., gGu, gGv, gGw]), for the function space functionSpaceType (e.g. "Lagrange" or "GradLagrange" for Lagrange basis functions or their gradient, in the u, v, w coordinates of the reference element; or "H1Legendre3" or "GradH1Legendre3" for 3rd order hierarchical H1 Legendre functions). numComponents returns the number C of components of a basis function. basisFunctions returns the value of the N basis functions at the evaluation points, i.e. [g1f1, g1f2, ..., g1fN, g2f1, ...] when C == 1 or [g1f1u, g1f1v, g1f1w, g1f2u, ..., g1fNw, g2f1u, ...] when C == 3. For basis functions that depend on the orientation of the elements, all values for the first orientation are returned first, followed by values for the second, etc. numOrientations returns the overall number of orientations. If wantedOrientations is not empty, only return the values for the desired orientation indices.

Return numComponents, basisFunctions, numOrientations.

source
gmsh.model.mesh.getBasisFunctionsOrientationForElementsFunction
gmsh.model.mesh.getBasisFunctionsOrientationForElements(elementType, functionSpaceType, tag = -1, task = 0, numTasks = 1)

Get the orientation index of the elements of type elementType in the entity of tag tag. The arguments have the same meaning as in getBasisFunctions. basisFunctionsOrientation is a vector giving for each element the orientation index in the values returned by getBasisFunctions. For Lagrange basis functions the call is superfluous as it will return a vector of zeros.

Return basisFunctionsOrientation.

source
gmsh.model.mesh.getNumberOfOrientationsFunction
gmsh.model.mesh.getNumberOfOrientations(elementType, functionSpaceType)

Get the number of possible orientations for elements of type elementType and function space named functionSpaceType.

Return an integer value.

source
gmsh.model.mesh.getEdgeNumberFunction
gmsh.model.mesh.getEdgeNumber(edgeNodes)

Get the global edge identifier edgeNum for an input list of node pairs, concatenated in the vector edgeNodes. Warning: this is an experimental feature and will probably change in a future release.

Return edgeNum.

source
gmsh.model.mesh.getLocalMultipliersForHcurl0Function
gmsh.model.mesh.getLocalMultipliersForHcurl0(elementType, tag = -1)

Get the local multipliers (to guarantee H(curl)-conformity) of the order 0 H(curl) basis functions. Warning: this is an experimental feature and will probably change in a future release.

Return localMultipliers.

source
gmsh.model.mesh.getKeysForElementsFunction
gmsh.model.mesh.getKeysForElements(elementType, functionSpaceType, tag = -1, returnCoord = true)

Generate the keys for the elements of type elementType in the entity of tag tag, for the functionSpaceType function space. Each key uniquely identifies a basis function in the function space. If returnCoord is set, the coord vector contains the x, y, z coordinates locating basis functions for sorting purposes. Warning: this is an experimental feature and will probably change in a future release.

Return keys, coord.

source
gmsh.model.mesh.getKeysForElementFunction
gmsh.model.mesh.getKeysForElement(elementTag, functionSpaceType, returnCoord = true)

Get the keys for a single element elementTag.

Return keys, coord.

source
gmsh.model.mesh.getNumberOfKeysForElementsFunction
gmsh.model.mesh.getNumberOfKeysForElements(elementType, functionSpaceType)

Get the number of keys by elements of type elementType for function space named functionSpaceType.

Return an integer value.

source
gmsh.model.mesh.getInformationForElementsFunction
gmsh.model.mesh.getInformationForElements(keys, elementType, functionSpaceType)

Get information about the keys. infoKeys returns information about the functions associated with the keys. infoKeys[0].first describes the type of function (0 for vertex function, 1 for edge function, 2 for face function and 3 for bubble function). infoKeys[0].second gives the order of the function associated with the key. Warning: this is an experimental feature and will probably change in a future release.

Return infoKeys.

source
gmsh.model.mesh.getBarycentersFunction
gmsh.model.mesh.getBarycenters(elementType, tag, fast, primary, task = 0, numTasks = 1)

Get the barycenters of all elements of type elementType classified on the entity of tag tag. If primary is set, only the primary nodes of the elements are taken into account for the barycenter calculation. If fast is set, the function returns the sum of the primary node coordinates (without normalizing by the number of nodes). If tag < 0, get the barycenters for all entities. If numTasks > 1, only compute and return the part of the data indexed by task.

Return barycenters.

source
gmsh.model.mesh.getElementEdgeNodesFunction
gmsh.model.mesh.getElementEdgeNodes(elementType, tag = -1, primary = false, task = 0, numTasks = 1)

Get the nodes on the edges of all elements of type elementType classified on the entity of tag tag. nodeTags contains the node tags of the edges for all the elements: [e1a1n1, e1a1n2, e1a2n1, ...]. Data is returned by element, with elements in the same order as in getElements and getElementsByType. If primary is set, only the primary (begin/end) nodes of the edges are returned. If tag < 0, get the edge nodes for all entities. If numTasks > 1, only compute and return the part of the data indexed by task.

Return nodeTags.

source
gmsh.model.mesh.getElementFaceNodesFunction
gmsh.model.mesh.getElementFaceNodes(elementType, faceType, tag = -1, primary = false, task = 0, numTasks = 1)

Get the nodes on the faces of type faceType (3 for triangular faces, 4 for quadrangular faces) of all elements of type elementType classified on the entity of tag tag. nodeTags contains the node tags of the faces for all elements: [e1f1n1, ..., e1f1nFaceType, e1f2n1, ...]. Data is returned by element, with elements in the same order as in getElements and getElementsByType. If primary is set, only the primary (corner) nodes of the faces are returned. If tag < 0, get the face nodes for all entities. If numTasks > 1, only compute and return the part of the data indexed by task.

Return nodeTags.

source
gmsh.model.mesh.getGhostElementsFunction
gmsh.model.mesh.getGhostElements(dim, tag)

Get the ghost elements elementTags and their associated partitions stored in the ghost entity of dimension dim and tag tag.

Return elementTags, partitions.

source
gmsh.model.mesh.setSizeFunction
gmsh.model.mesh.setSize(dimTags, size)

Set a mesh size constraint on the model entities dimTags. Currently only entities of dimension 0 (points) are handled.

source
gmsh.model.mesh.setSizeAtParametricPointsFunction
gmsh.model.mesh.setSizeAtParametricPoints(dim, tag, parametricCoord, sizes)

Set mesh size constraints at the given parametric points parametricCoord on the model entity of dimension dim and tag tag. Currently only entities of dimension 1 (lines) are handled.

source
gmsh.model.mesh.setSizeCallbackFunction
gmsh.model.mesh.setSizeCallback(callback)

Set a global mesh size callback. The callback should take 5 arguments (dim, tag, x, y and z) and return the value of the mesh size at coordinates (x, y, z).

source
gmsh.model.mesh.setTransfiniteCurveFunction
gmsh.model.mesh.setTransfiniteCurve(tag, numNodes, meshType = "Progression", coef = 1.)

Set a transfinite meshing constraint on the curve tag, with numNodes nodes distributed according to meshType and coef. Currently supported types are "Progression" (geometrical progression with power coef), "Bump" (refinement toward both extremities of the curve) and "Beta" (beta law).

source
gmsh.model.mesh.setTransfiniteSurfaceFunction
gmsh.model.mesh.setTransfiniteSurface(tag, arrangement = "Left", cornerTags = Cint[])

Set a transfinite meshing constraint on the surface tag. arrangement describes the arrangement of the triangles when the surface is not flagged as recombined: currently supported values are "Left", "Right", "AlternateLeft" and "AlternateRight". cornerTags can be used to specify the (3 or 4) corners of the transfinite interpolation explicitly; specifying the corners explicitly is mandatory if the surface has more that 3 or 4 points on its boundary.

source
gmsh.model.mesh.setTransfiniteVolumeFunction
gmsh.model.mesh.setTransfiniteVolume(tag, cornerTags = Cint[])

Set a transfinite meshing constraint on the surface tag. cornerTags can be used to specify the (6 or 8) corners of the transfinite interpolation explicitly.

source
gmsh.model.mesh.setTransfiniteAutomaticFunction
gmsh.model.mesh.setTransfiniteAutomatic(dimTags = Tuple{Cint,Cint}[], cornerAngle = 2.35, recombine = true)

Set transfinite meshing constraints on the model entities in dimTag. Transfinite meshing constraints are added to the curves of the quadrangular surfaces and to the faces of 6-sided volumes. Quadragular faces with a corner angle superior to cornerAngle (in radians) are ignored. The number of points is automatically determined from the sizing constraints. If dimTag is empty, the constraints are applied to all entities in the model. If recombine is true, the recombine flag is automatically set on the transfinite surfaces.

source
gmsh.model.mesh.setRecombineFunction
gmsh.model.mesh.setRecombine(dim, tag)

Set a recombination meshing constraint on the model entity of dimension dim and tag tag. Currently only entities of dimension 2 (to recombine triangles into quadrangles) are supported.

source
gmsh.model.mesh.setSmoothingFunction
gmsh.model.mesh.setSmoothing(dim, tag, val)

Set a smoothing meshing constraint on the model entity of dimension dim and tag tag. val iterations of a Laplace smoother are applied.

source
gmsh.model.mesh.setReverseFunction
gmsh.model.mesh.setReverse(dim, tag, val = true)

Set a reverse meshing constraint on the model entity of dimension dim and tag tag. If val is true, the mesh orientation will be reversed with respect to the natural mesh orientation (i.e. the orientation consistent with the orientation of the geometry). If val is false, the mesh is left as-is.

source
gmsh.model.mesh.setAlgorithmFunction
gmsh.model.mesh.setAlgorithm(dim, tag, val)

Set the meshing algorithm on the model entity of dimension dim and tag tag. Currently only supported for dim == 2.

source
gmsh.model.mesh.setSizeFromBoundaryFunction
gmsh.model.mesh.setSizeFromBoundary(dim, tag, val)

Force the mesh size to be extended from the boundary, or not, for the model entity of dimension dim and tag tag. Currently only supported for dim ==

source
gmsh.model.mesh.setCompoundFunction
gmsh.model.mesh.setCompound(dim, tags)

Set a compound meshing constraint on the model entities of dimension dim and tags tags. During meshing, compound entities are treated as a single discrete entity, which is automatically reparametrized.

source
gmsh.model.mesh.setOutwardOrientationFunction
gmsh.model.mesh.setOutwardOrientation(tag)

Set meshing constraints on the bounding surfaces of the volume of tag tag so that all surfaces are oriented with outward pointing normals. Currently only available with the OpenCASCADE kernel, as it relies on the STL triangulation.

source
gmsh.model.mesh.embedFunction
gmsh.model.mesh.embed(dim, tags, inDim, inTag)

Embed the model entities of dimension dim and tags tags in the (inDim, inTag) model entity. The dimension dim can 0, 1 or 2 and must be strictly smaller than inDim, which must be either 2 or 3. The embedded entities should not be part of the boundary of the entity inTag, whose mesh will conform to the mesh of the embedded entities.

source
gmsh.model.mesh.removeEmbeddedFunction
gmsh.model.mesh.removeEmbedded(dimTags, dim = -1)

Remove embedded entities from the model entities dimTags. if dim is >= 0, only remove embedded entities of the given dimension (e.g. embedded points if dim == 0).

source
gmsh.model.mesh.reorderElementsFunction
gmsh.model.mesh.reorderElements(elementType, tag, ordering)

Reorder the elements of type elementType classified on the entity of tag tag according to ordering.

source
gmsh.model.mesh.setPeriodicFunction
gmsh.model.mesh.setPeriodic(dim, tags, tagsMaster, affineTransform)

Set the meshes of the entities of dimension dim and tag tags as periodic copies of the meshes of entities tagsMaster, using the affine transformation specified in affineTransformation (16 entries of a 4x4 matrix, by row). If used after meshing, generate the periodic node correspondence information assuming the meshes of entities tags effectively match the meshes of entities tagsMaster (useful for structured and extruded meshes). Currently only available for @code{dim} == 1 and @code{dim} == 2.

source
gmsh.model.mesh.getPeriodicNodesFunction
gmsh.model.mesh.getPeriodicNodes(dim, tag, includeHighOrderNodes = false)

Get the master entity tagMaster, the node tags nodeTags and their corresponding master node tags nodeTagsMaster, and the affine transform affineTransform for the entity of dimension dim and tag tag. If includeHighOrderNodes is set, include high-order nodes in the returned data.

Return tagMaster, nodeTags, nodeTagsMaster, affineTransform.

source
gmsh.model.mesh.splitQuadranglesFunction
gmsh.model.mesh.splitQuadrangles(quality = 1., tag = -1)

Split (into two triangles) all quadrangles in surface tag whose quality is lower than quality. If tag < 0, split quadrangles in all surfaces.

source
gmsh.model.mesh.classifySurfacesFunction
gmsh.model.mesh.classifySurfaces(angle, boundary = true, forReparametrization = false, curveAngle = pi, exportDiscrete = true)

Classify ("color") the surface mesh based on the angle threshold angle (in radians), and create new discrete surfaces, curves and points accordingly. If boundary is set, also create discrete curves on the boundary if the surface is open. If forReparametrization is set, create edges and surfaces that can be reparametrized using a single map. If curveAngle is less than Pi, also force curves to be split according to curveAngle. If exportDiscrete is set, clear any built-in CAD kernel entities and export the discrete entities in the built- in CAD kernel.

source
gmsh.model.mesh.createGeometryFunction
gmsh.model.mesh.createGeometry(dimTags = Tuple{Cint,Cint}[])

Create a geometry for the discrete entities dimTags (represented solely by a mesh, without an underlying CAD description), i.e. create a parametrization for discrete curves and surfaces, assuming that each can be parametrized with a single map. If dimTags is empty, create a geometry for all the discrete entities.

source
gmsh.model.mesh.createTopologyFunction
gmsh.model.mesh.createTopology(makeSimplyConnected = true, exportDiscrete = true)

Create a boundary representation from the mesh if the model does not have one (e.g. when imported from mesh file formats with no BRep representation of the underlying model). If makeSimplyConnected is set, enforce simply connected discrete surfaces and volumes. If exportDiscrete is set, clear any built-in CAD kernel entities and export the discrete entities in the built-in CAD kernel.

source
gmsh.model.mesh.computeHomologyFunction
gmsh.model.mesh.computeHomology(domainTags = Cint[], subdomainTags = Cint[], dims = Cint[])

Compute a basis representation for homology spaces after a mesh has been generated. The computation domain is given in a list of physical group tags domainTags; if empty, the whole mesh is the domain. The computation subdomain for relative homology computation is given in a list of physical group tags subdomainTags; if empty, absolute homology is computed. The dimensions homology bases to be computed are given in the list dim; if empty, all bases are computed. Resulting basis representation chains are stored as physical groups in the mesh.

source
gmsh.model.mesh.computeCohomologyFunction
gmsh.model.mesh.computeCohomology(domainTags = Cint[], subdomainTags = Cint[], dims = Cint[])

Compute a basis representation for cohomology spaces after a mesh has been generated. The computation domain is given in a list of physical group tags domainTags; if empty, the whole mesh is the domain. The computation subdomain for relative cohomology computation is given in a list of physical group tags subdomainTags; if empty, absolute cohomology is computed. The dimensions homology bases to be computed are given in the list dim; if empty, all bases are computed. Resulting basis representation cochains are stored as physical groups in the mesh.

source
gmsh.model.mesh.computeCrossFieldFunction
gmsh.model.mesh.computeCrossField()

Compute a cross field for the current mesh. The function creates 3 views: the H function, the Theta function and cross directions. Return the tags of the views.

Return viewTags.

source
gmsh.model.mesh.triangulateFunction
gmsh.model.mesh.triangulate(coord)

Triangulate the points given in the coord vector as pairs of u, v coordinates, and return the node tags (with numbering starting at 1) of the resulting triangles in tri.

Return tri.

source
gmsh.model.mesh.tetrahedralizeFunction
gmsh.model.mesh.tetrahedralize(coord)

Tetrahedralize the points given in the coord vector as triplets of x, y, z coordinates, and return the node tags (with numbering starting at 1) of the resulting tetrahedra in tetra.

Return tetra.

source

Mesh Size

gmsh.model.mesh.field.addFunction
gmsh.model.mesh.field.add(fieldType, tag = -1)

Add a new mesh size field of type fieldType. If tag is positive, assign the tag explicitly; otherwise a new tag is assigned automatically. Return the field tag.

Return an integer value.

source