Manipulation of meshes in 2-D

    

We detail here the attributes and methods related to the class Mesh<Dimension2>, they are sorted by order of appearance for each class. We call "referenced edges" all the edges which are explicitly provided in the mesh, usually these edges are located on the boundary of the domain or at the interface between two domains. Boundary conditions will be applied on "referenced edges", so if you want to specify a boundary condition to a curve, you need to have it referenced. Each referenced edge possess a reference number, which is an identification number of the boundary to which the edge belongs. Two elementary entities are detailed : Edge<Dimension2> and Face<Dimension2> .

Methods for Edge<Dimension2>

Init sets the two extremity numbers of the edge
SetReference sets the reference number of the edge
GetReference returns the reference number of the edge
SetGeometryReference sets the geometry number of the edge
GetGeometryReference returns the geometry number of the edge
numVertex returns vertex number of an extremity of the edge
GetNbVertices returns 2
ClearConnectivity sets element numbers to -1
GetNbElements returns the number of elements sharing the edge
numElement returns element number of an element adjacent to the edge
AddElement adds an element sharing this edge
SetElement changes an element sharing this edge

Methods for Face<Dimension2>

IsTriangular returns true if the face is a triangle
IsQuadrangular returns true if the face is a quadrangle
GetHybridType returns 0 for a triangle and 1 for a quadrangle
GetOrientationEdge returns the difference of orientation between local and global edge
GetOrientationBoundary returns the difference of orientation between local and global edge
numVertex returns vertex number of a vertex of the face
numEdge returns edge number of an edge of the face
GetReference returns reference number associated with the face
GetNbVertices returns the number of vertices of the face
GetNbEdges returns the number of edges of the face
SetVertex sets vertex number of a vertex of the face
SetReference sets reference number associated with the face
InitTriangular sets all three vertex numbers of the face
InitQuadrangular sets all four vertex numbers of the face
Init sets all vertex numbers of the face
SetEdge sets an edge of the face
GetPositionBoundary returns local position of an edge in the face
ClearConnectivity removes information about edges
ClearEdges removes information about edges
IsPML returns true if the element is PML
SetPML sets the element as a PML element
UnsetPML sets the element as a non-PML element
GetMemorySize returns the size used by the object in bytes
GetTypePML returns the type of the PML associated with the face
GetNumberPML returns the number of the PML region associated with the face
IsCurved returns true if the element is curved
SetCurved sets the element as a curved element
UnsetCurved sets the element as a non-curved element
GetNbBoundary returns the number of edges of the face
numBoundary returns edge number

Below, we detail public attributes of a mesh :

Public attributes (class Mesh<Dimension2>)

print_level level of messages to be printed (-1 for silent mode, 10 for verbose mode)
GlobElementNumber_Subdomain Global element numbers
GlobFaceNumber_Subdomain Global face numbers
GlobEdgeNumber_Subdomain Global edge numbers
GlobVertexNumber_Subdomain Global vertex numbers
Glob_invSqrtJacobian inverse of square root of jacobian on quadrature and nodal points
Glob_invSqrtJacobianNeighbor inverse of square root of jacobian on quadrature points of neighbor elements
Glob_GradJacobian gradient of jacobian on quadrature points

Public attributes (class MeshNumbering<Dimension2>)

number_map numbering scheme used to number the mesh
OffsetDofFaceNumber Offsets for dof numbers on the interior of each face (triangle/quadrangle)
OffsetDofEdgeNumber Offsets for dof numbers on the interior of each edge
OffsetDofElementNumber Offsets for dof numbers on the interior of each element (triangle/quadrangle)
OffsetDofVertexNumber Offsets for dof numbers on vertices
OffsetQuadElementNumber Offsets for quadrature points on the boundaries of elements
GlobDofNumber_Subdomain Global dof numbers
GlobDofPML_Subdomain Global dof numbers for PML layers
treat_periodic_condition_during_number if true, periodic dofs are searched when the mesh is numbered
compute_dof_pml if true, PML dofs are computed
drop_interface_pml_dof if true, dofs at the interface between the physical and PMLs are not PML dofs

Then we detail all the methods related to meshes, we sort the method by the class where it is defined.

Methods specific to the class Mesh<Dimension2> (inherited from Mesh_Base<Dimension2> and MeshBoundaries<Dimension2>)

Boundary gives access to an edge of the mesh
GetNbBoundary returns the number of edges in the mesh
ReallocateBoundariesRef changes the size of the array containing the referenced edges (previous edges are lost)
ResizeBoundariesRef changes the size of the array containing the referenced edges (previous edges are kept)
BoundaryRef gives access to a referenced edge of the mesh
GetNbBoundaryRef returns the number of referenced edges in the mesh
GetNbEdgesRef returns the number of referenced edges in the mesh
GetTriangleReferenceElement returns triangular nodal finite element used for definition of curved triangles
GetQuadrilateralReferenceElement returns quadrilateral nodal finite element used for definition of curved quadrilaterals
GetReferenceElementSurface returns the 1-D nodal finite element
GetReferenceElementVolume returns the list of nodal finite elements used for each element
GetLocalNodalNumber returns local numbers of nodes of an edge
HybridMesh returns true if the mesh contains quadrilaterals and triangles
IsOnlyTriangular returns true if the mesh contains only triangles
IsOnlyQuadrilateral returns true if the mesh contains only quadrilaterals
GetMemorySize returns the memory used to store the object in bytes
GetNbTriangles returns the number of triangles in the mesh
GetNbQuadrangles returns the number of quadrilaterals in the mesh
SetGeometryOrder defines the order of approximation for the geometry (curved mesh)
GetUpperBoundNumberOfEdges returns an upper bound of the number of edges in the mesh
SetInputData modification of the object with an entry of the data file.
GetNodesCurvedMesh computes all the internal nodes of a high order mesh, and the associated numbering
Read reads a mesh from a file
Write writes a mesh on a file
ReorientElements permutes local numbering of elements in order to have only positive jacobians.
FindConnectivity finds all the edges of the mesh and links them with referenced edges and elements
RemoveDuplicateVertices sorts vertices and removes duplicate vertices
ForceCoherenceMesh removes unused vertices and referenced edges
ClearZeroBoundaryRef removes boundary edges whose reference is equal to 0
ApplyLocalRefinement performs a local refinement around a vertex of the mesh
AddVertexConformal adds a new vertex in the mesh (the mesh is kept conformal)
SymmetrizeMesh appends to the current mesh the image obtained by a symmetry
PeriodizeMeshTeta generates a cyclic mesh from the section
SplitIntoQuadrilaterals splits each triangle into 3 quadrileterals, and each quadrilateral into 4 quadrilaterals
SplitIntoTriangles splits each quadrilateral into two triangles
GetMeshSubdivision computes all the nodes of the high-order mesh and numbering for each element
SubdivideMesh splits each element into small elements (triangles into small triangles, quadrangles into small quadrangles)
SymmetrizeMeshToGetPositiveX If the mesh has only negative x-coordinate, it is symmetrized to have positive x-coordinate
CreateRegularMesh creation of a mesh of a rectangle with quadrangles or triangles.
ExtrudeCoordinate adds layers of elements along a coordinate (x or y) when the outer boundary is a rectangle
ConstructMesh constructs the mesh with the instructions given in the data file
RearrangeElements changes the ordering of elements such that triangles are placed first
CreateSubmesh creates a mesh from a set of elements of the current mesh
AppendMesh appends a mesh to the current mesh
GetEdgesOnZaxis finds all the edges on axis x = 0
GetTetaMinMax finds minimal and maximal value of teta (polar coordinates) in the mesh

Methods from MeshNumbering<Dimension2>

CopyInputData input parameters are copied from another mesh numbering
ReallocateNeighborElements allocates arrays to store the orders fo neighbor elements (parallel)
FinalizeNeighborElements finalizes arrays to store the orders fo neighbor elements (parallel)
GetMemorySize returns the size used by the object in bytes
SetInputData modification of the object with an entry of the data file.
ClearDofs removes dof numbers from the mesh
Clear clears the numbering
NumberMesh clears previous numbers and renumbers the mesh
Element gives access to the numbering of an element
GetNbNeighborElt returns the number of neighboring elements
GetOrderNeighborElement returns the order of a neighboring element
GetTypeNeighborElement returns the type of a neighboring element
GetEdgeNeighborElement returns the edge number associated with a neighboring element
GetLocalEdgeNumberNeighborElement returns the neighboring element number from the edge number
SetOrderNeighborElement sets the order of a neighboring element
SetTypeNeighborElement sets the type of a neighboring element
SetEdgeNeighborElement sets the edge number associated with a neighboring element
GetPeriodicityTypeForDof returns the type of periodicity (in theta, x, y, z) for a periodic dof
SetPeriodicityTypeForDof sets the type of periodicity (in theta, x, y, z) for a periodic dof
GetPeriodicityTypeForBoundary returns the type of periodicity (in theta, x, y, z) for a periodic boundary
SetPeriodicityTypeForBoundary sets the type of periodicity (in theta, x, y, z) for a periodic boundary
GetNbPeriodicDof returns the number of dofs on periodic boundaries
GetPeriodicDof returns dof number obtained by a translation of an original dof
SetPeriodicDof sets dof number obtained by a translation of an original dof
GetOriginalPeriodicDof returns original dof number
SetOriginalPeriodicDof sets original dof number
GetPeriodicBoundary returns the number of the periodic boundary
SetPeriodicBoundary sets the number of the periodic boundary
SetSameNumberPeriodicDofs same dof numbers are used (continuous elements)
SetFormulationForPeriodicCondition informs how periodic condition is enforced
GetFormulationForPeriodicCondition returns how periodic condition is enforced
UseSameNumberForPeriodicDofs returns true if same dof numbers are used for periodic condition
GetOrder returns the order(s) of approximation associated with the mesh
GetNbDof return the number of degrees of freedom in the mesh
SetNbDof sets the number of degrees of freedom in the mesh
GetNbLocalDof returns the number of degrees of freedom of an element
GetNbDofPML returns the number of degrees of freedom in PMLs
GetDofPML returns the pml dof number from the dof number
ReallocateDofPML sets the number of degrees of freedom in PMLs
SetDofPML sets a pml dof number
SetOrder sets the order of approximation associated with the mesh
IsOrderVariable returns true if the order is variable (different orders used in the mesh)
SetVariableOrder sets the strategy used to affect orders
SetMeshSizeVariableOrder sets the mesh sizes for each order in the case of variable order
SetCoefficientVariableOrder sets a coefficient (for each domain) for the mesh sizes in the case of variable order
GetOrderElement returns the order of approximation associated with an element of the mesh
GetOrderEdge returns the order of approximation associated with an edge of the mesh
GetOrderFace returns the order of approximation associated with a face of the mesh
GetOrderInside returns the order of approximation associated with the interior of an element
SetOrderEdge sets the order of approximation associated with an edge of the mesh
SetOrderFace sets the order of approximation associated with a face of the mesh
SetOrderInside sets the order of approximation associated with the interior of an element of the mesh
GetNbPointsQuadratureBoundary returns the number of quadrature points on a boundary of the mesh
InitPeriodicBoundary initializes periodic boundaries
GetTranslationPeriodicDof returns translation vector between two dofs
SetTranslationPeriodicDof sets translation vector between two dofs
GetTranslationPeriodicBoundary returns translation vector between two boundaries
ReallocateElements changes the number of elements stored in the numbering.
ReallocatePeriodicDof changes the number of periodic dofs
ComputeVariableOrder computes order for each edge and element
GetLocalVertexBoundary returns the vertex number of an edge of the reference element
GetOrderQuadrature returns the order of quadrature formula associated with an edge of the mesh
SetOrderQuadrature sets the order of quadrature associated with an edge of the mesh
GetReferenceQuadrature returns quadrature points used for the unit interval
SetOrderElement sets the order of approximation associated with an element of the mesh
GetRotationFace returns the difference of orientation between two edges
GetOppositeOrientationFace returns the opposite orientation between two edges
GetNbDofElement returns the number of dofs for an element by using informations of number_map
GetBoundaryRotation computation of difference of orientation between two elements sharing a common edge.
ReconstructOffsetDofs recomputes offsets for dofs on edges, faces, and elements
GetNbDofBoundary return the number of dofs for a boundary of the mesh
GetDofBoundary retrieves the numbers of dofs associated with a boundary
ConstructQuadrilateralNumbering returns numbering scheme used for nodes in a quadrilateral
ConstructTriangularNumbering returns numbering scheme used for nodes in a triangle
TreatPeriodicCondition finds corresponding dofs for periodic/cyclic condition

Methods from MeshBoundaries<Dimension2>

GetNbReferences returns the number of references handled in the mesh
GetCurveType returns the type of the curve (circle, spline, etc) associated with a reference
GetBoundaryCondition returns the boundary condition (Dirichlet, Neumann, Robin, etc) associated with a reference
GetBodyNumber returns the body number associated with a reference
GetCurveParameter returns the parameters of the curve (for a circle, radius and center) associated with a reference
GetCurveFunctions1D returns the basis functions associated with 1-D nodal points
GetNbPmlAreas returns the number of PML areas
ReallocatePmlAreas reallocates the PML areas
GetPmlArea gives access to a PML area
SetBoundaryCondition sets the boundary condition (Dirichlet, Neumann, Robin, etc) associated with a reference
SetCurveType sets the type of the curve (circle, spline, etc) associated with a reference
SetBodyNumber sets the body number associated with a reference
SetCurveParameter sets the parameters of the curve (for a circle, radius and center) associated with a reference
CopyCurves curve parameters are copied from another mesh
CopyInputData input parameters are copied from another mesh
GetNewReference returns an unused reference so that you can create referenced edges with that reference and specify boundary condition or curved boundary.
ResizeNbReferences changes the number of references handled in the mesh
SetNewReference sets boundary condition and curve type associated with a reference
GetNbPeriodicReferences returns the number of periodic boundaries
GetPeriodicReference returns the reference of periodic boundaries
AddPeriodicCondition adds periodic boundaries
GetPeriodicityTypeReference returns type of periodicity (in theta, x, y, z) for a couple of references
GetPeriodicAlpha returns angle of cyclic section (2pi/N)
SetPeriodicAlpha sets angle of cyclic section (2pi/N)
ClearPeriodicCondition clears informations about periodic conditions
GetSafetyCoef returns safety coefficient to use for a curved element
AreCurveEqual returns true if two curves are equal
GetGeometryOrder returns the order used to approximate geometry
GetInterpolate1D returns value of 1-D shape function used for the approximation of geometry
GetNodalPointInsideEdge returns intern points used on the edge [0, 1] for the geometry (should be Gauss-Lobatto points)
FindEdgeRefWithExtremities returns referenced edge number, knowing the two extremities of the edge
FindEdgeWithExtremities returns edge number, knowing the two extremities of the edge
ClearCurves clears informations about curves
IsBoundaryCurved returns true if the edge i is curved
GetPointInsideEdge returns intermediary point on a curved edge
SetPointInsideEdge modifies intermediary point on a curved edge
FindFollowingVertex returns the next point on a curve, knowing two point on this curve (so that you can loop on the curve)
GetSurfaceFiniteElement returns the 1-D finite element used for curves
IsPointOnBoundaryRef returns true if a point is on a referenced edge
GetAllReferencesAroundPoint returns all the references of edges sharing a given point
GetPointWithOtherReference returns a vertex near another given vertex, but belonging to a different reference
RemoveReference removes referenced edges with a given reference
SortPeriodicBoundaries modifies vertex numbers in order to have the same orientations between periodic edges
SortBoundariesRef modifies referenced edges so that for each edge, the first extremity has a lower number than second extremity
ConstructCrossReferenceBoundaryRef links referenced edges with global edges
AddBoundaryEdges finds edges on the boundary which are not already in referenced edges, and adds them
ProjectPointsOnCurves finds parameters of curves, projects points on them, and computes position of internal points on curved boundaries
RedistributeReferences redistributes references numbers such that each closed curve has a different number
CheckCurveParameter checks that parameters of a given curve seem correct
ReadCurveType sets the type of curve with its name (string coming from the data file)
ProjectToCurve projects a point on a curve
GetPointsOnCurve computes a point on a parametrized curve for a given value of the parameter
GetAllPointsOnReference return all the points (and normales) belonging to referenced edges with the same reference
GetDistantPointsOnReference return some points (and normales) belonging to referenced edges with the same reference
FindParametersPlane finds coefficients of a straight line
FindParametersSpline finds coefficients of a spline
FindParametersCircle finds coefficients of a circle
FindParametersEllipse finds coefficients of an ellipse
FindParametersPeanut finds coefficients of a peanut
FindParametersCurve finds parameters of all the curves (for a circle, radius and center)
GetCenterRadialPML returns the center of the radial PMLs (if present)
GetBoundaryMesh extracts some referenced boundaries of the current mesh and creates a new mesh with that
GenerateSurfaceOfRevolution Generates a 3-D surfacic mesh (only triangles) by revolution of a 2-D curve around z-axis

Methods of PmlRegionParameter

GetThicknessPML returns the thickness of PML
GetOriginRadialPML returns the center of radial PML
GetRadiusPML returns the radius of radial PML
SetThicknessPML sets the thickness of PML
GetTypeAdditionPML returns the type of PML to add to the domain
GetNbLayersPML returns the number of cells with respect to PML direction
SetRadiusPML sets the radius of radial PML
SetOriginRadialPML sets the center of radial PML
SetAdditionPML sets PML to add to the domain
SetParameters sets PML with a line of the data file
AddPML adds PML elements to a mesh
GetAverageLengthOnBoundary returns the mesh size on the outer boundary
AddPMLElements adds PML elements on the outer boundary of the mesh
GetNbParameters returns the number of parameters associated with the PML
FillParameters fills an array with parameters associated with the PML
SetRegion sets the pml with an input array

Methods from Mesh_Base<Dimension2>

GetNbElt returns the number of elements in the mesh
GetNbVertices returns the number of vertices in the mesh
GetNbEdges returns the number of edges in the mesh
GetOriginalNeighborReference returns original references of faces shared by two processors
SetOriginalNeighborReference sets original references of faces shared by two processors
GetPathName returns path name used in method ConstructMesh
SetPathName changes path name used in method ConstructMesh
IsElementAffine returns true if the considered element is obtained from the reference element by an affine transformation
ReallocateElements changes the number of elements stored in the mesh. Previous elements are not kept.
ResizeElements changes the number of elements stored in the mesh. Previous elements are kept.
ClearElements removes all the elements of the mesh
Element gives access to an element of the mesh
GetTypeElement returns the hybrid nature of an element (0 for triangle, 1 for quadrangle)
GetEdge gives access to an edge of the mesh
ReallocateEdges changes the number of edges stored in the mesh. Previous edges are not kept.
ResizeEdges changes the number of edges stored in the mesh. Previous edges are kept.
ReallocateVertices changes the number of vertices stored in the mesh. Previous vertices are lost.
ResizeVertices changes the number of vertices stored in the mesh. Previous vertices are kept.
ClearVertices removes all the vertices
Vertex gives access to a vertex of the mesh
GetXmin returns the minimum of x-coordinate among vertices of the mesh
GetYmin returns the minimum of y-coordinate among vertices of the mesh
GetZmin returns the minimum of z-coordinate among vertices of the mesh
GetXmax returns the maximum of x-coordinate among vertices of the mesh
GetYmax returns the maximum of y-coordinate among vertices of the mesh
GetZmax returns the maximum of z-coordinate among vertices of the mesh
FindVertexNumber returns the number of a vertex from its coordinates
GetReferenceElement returns the finite element for the i-th element
GetVerticesElement retrieves the vertices coordinates of an element
GetDofsElement retrieves the positions of degrees of freedom of an element
SetInputData modification of object with an entry of the data file
PointsOnSameEdge returns true if two points belong to the same edge
GetNeighboringElementsAroundVertex finds all the elements around a given point
GetNeighboringEdgesAroundVertex finds all the edges around a given point
SetVerticesToBeAdded sets the points that will be added to the mesh in ConstructMesh
SetVerticesToBeRefined sets the points that will need a local refinement in ConstructMesh
CreateReferenceVertices attributes a reference to each vertex of the mesh
GetMeshSize returns an evaluation of the mesh size
ClearConnectivity removes all the edges of the mesh (but not referenced edges)
ClearInput clears input parameters used in ConstructMesh
Clear clears the mesh
AddOverlapRegion appends layers of elements to a region of the mesh
SplitIntoBoxes splits a mesh into several subdomains by using a regular cartesian grid
SplitMetis splits a mesh into several subdomains by using Metis
SplitScotch splits a mesh into several subdomains by using Scotch
SplitConcentric splits a mesh into several subdomains by using a basic concentric algorithm
SplitLayered splits a mesh into several subdomains by using references
PartMesh splits a mesh into several subdomains by providing directly epart
UpdatePeriodicReferenceSplit updates the splitting of the mesh by taking into account periodic condition
WriteOrder writes a mesh in a file, while replacing physical references of element by orders
FjElemNodal computes nodal points for an element
DFjElemNodal computes jacobian matrix on nodal points for an element
Fj computes images of points by transformation Fi
DFj computes jacobian matrix for an element
GetNbPointsNodalElt returns the number of nodal points of an element
GetNbPointsNodalBoundary returns the number of nodal points on an edge of a face
GetNodalNumber returns the number of a nodal point from its number on an edge
GetNormale computes the normale from the jacobian matrix
GetDistanceToBoundary returns the distance to the boundary of the reference element
ProjectPointOnBoundary projects a point to the boundary of the reference element
GetBoundingBox computes a rectangular box bounding an element
ComputeValuesPhiNodalRef computation of the shape functions on a given point
ComputeNodalGradient computation of the gradient from values on nodal points of an element

Functions related to meshes in 2-D

ExtrudePMLLayer adds elements for PML in the case of slanted/circular PML
TranslateMesh Applies a translation to the mesh
RotateMesh Applies a rotation to the mesh
ScaleMesh Each coordinate (x, y) is scaled differently
WriteElementMesh Writes a single element of the mesh in a file
SendMeshToProcessor Sends a mesh to another processor
RecvMeshFromProcessor Receives a mesh from another processor

Surfacic meshes

Surface meshes are meshes of a boundary of an initial mesh. It is usually constructed by calling the method GetBoundaryMesh. It consists essentially of a mesh supplemented by the arrays ListeBoundaries, ListeVertices, NumElement and NumLocalBoundary. Below we show an exemple of use

  // the "volume" mesh
  Mesh<Dimension2> mesh_vol;
  mesh_vol.Read("init.mesh");

  // the "surface" mesh
  SurfacicMesh<Dimension2> mesh_surf;
  // we want to extract the boundaries of reference 2 and 3 for instance :
  IVect ref_cond(mesh_vol.GetNbReferences()+1); ref_cond.Zero();
  ref_cond(2) = 1; ref_cond(3) = 1;

  // we call GetBoundaryMesh
  mesh_vol.GetBoundaryMesh(1, mesh_surf, ref_cond);

  // in mesh_surf, you have the list of extracted edges : mesh_surf.ListeBoundaries
  DISP(mesh_surf.ListeBoundaries);
  // the list of extracted vertices : mesh_surf.ListeVertices
  DISP(mesh_surf.ListeVertices);

  // for each extracted edge, you know the element number : mesh_surf.NumElement
  DISP(mesh_surf.NumElement);

  // and the local position of the extracted edge in this element
  DISP(mesh_surf.NumLocalBoundary);

print_level

This attribute stores the level of messages to be printed. You can set it to -1 for silent mode (nothing will be displayed), to 10 for a very verbose mode, 2 for an intermediate mode.

Example :

  Mesh<Dimension2> mesh;
  mesh.print_level = 3;

Location :

MeshBase.hxx

GlobElementNumber_Subdomain, GlobFaceNumber_Subdomain, GlobEdgeNumber_Subdomain, GlobVertexNumber_Subdomain

These attributes are used in Montjoie to store the global element number, face number, edge number and vertex number. They are relevant only for parallel computation where each processor stores his part of the mesh, you can access to global numbers with these arrays.

Example :

  Mesh<Dimension2> mesh;

  // global number of element 0 ?
  int num_glob = mesh.GlobElementNumber_Subdomain(0);

Location :

MeshBase.hxx

Glob_invSqrtJacobian, Glob_invSqrtJacobianNeighbor, Glob_GradJacobian

These attributes are used in Montjoie to store the inverse of square root of Ji (determinant of DFi). It is used for the Warburton's trick (useful for curved tetrahedra).

Location :

MeshBase.hxx

number_map

This attribute stores the numbering scheme to use, i.e. number of dofs per vertex, edge, triangle, tetrahedron, how the numbers are transformed after rotation of a face. If you want more details, you can look at the methods of class NumberMap.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
// constant order on all the mesh
mesh_num.SetOrder(order);

// number of dofs per vertex, edge, etc
mesh_num.number_map.SetNbDofVertex(order, 1);
mesh_num.number_map.SetNbDofEdge(order, order-1);
mesh_num.number_map.SetNbDofQuadrangle(order, (order-1)*(order-1));
mesh_num.number_map.SetNbDofTriangle(order, (order-1)*(order-2)/2);

// then you can number the mesh
mesh_num.NumberMesh();

// and display global numbers
for (int i = 0; i < mesh.GetNbElt(); i++)
  for (int j = 0; j < mesh_num.Element(i).GetNbDof(); j++)
    cout << "Global dof number of local dof " << j
         << " and element  " << i << " : " << mesh_num.Element(i).GetNumberDof(j) << endl;

Location :

NumberMesh.hxx

OffsetDofVertexNumber, OffsetDofEdgeNumber, OffsetDofFaceNumber, OffsetDofElementNumber, OffsetQuadElementNumber

The attribute OffsetDofVertexNumber (continuous elements only) stores offsets used to number dofs on vertices. The attribute OffsetDofEdgeNumber (continuous elements only) stores offset used to number dofs on edges. For continuous elements, the attribute OffsetDofElementNumber stores offset used to number interior dofs on elements. Fod discontinuous elements, this attribute stores offset used to number all the dofs (that are only associated with elements. The attribute OffsetDofFaceNumber (HDG formulation only) stores offset used to number boundary unknown λ on edges of the mesh. Finally the attribute OffsetQuadElementNumber stores the offsets for the quadrature points on edges of each element.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// we consider a 2-D mesh numbering (with continuous elements)
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// if you want to know the first dof associated with a given vertex :
int num_vertex = 15;
int first_dof = mesh_num.OffsetDofVertexNumber(num_vertex);

// if you want to know the first dof associated with a given edge :
int num_edge = 21;
int first_dof_edge = mesh_num.OffsetDofVertexNumber(num_edge);

// if you want to know the first dof associated with a given element :
int num_elem = 3;
int first_dof_elt = mesh_num.OffsetDofElementNumber(num_elem);

// offset for quadrature points
int offset_quad = mesh_num.OffsetQuadElementNumber(num_elem);

Location :

NumberMesh.hxx

GlobDofNumber_Subdomain

The attribute GlobDofNumber_Subdomain stores the global dof number (for parallel computation) of a local dof. The attribute GlobDofPML_Subdomain stores the global pml dof numbers of local pml dofs

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// we consider a 2-D mesh numbering (with continuous elements)
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// global dof of a local dof (assuming that the mesh has been split and distributed over processors)
int num_loc = 42;
int glob_dof = mesh_num.GlobDofNumber_Subdomain(num_loc);

// and for PMLs :
int num_loc_pml = 12;
int glob_dof_pml = mesh_num.GlobDofPML_Subdomain(num_loc_pml);

Location :

NumberMesh.hxx

treat_periodic_condition_during_number

If the attribute treat_periodic_condition_during_number is true, the periodic conditions will be treated when the mesh is numbered. Otherwise the periodic conditions are not treated. By default this attribute is set to true.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// we consider a 2-D mesh numbering (with continuous elements)
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.treat_periodic_condition_during_number = false; // if you do not want to find periodic dofs
mesh_num.NumberMesh();

Location :

NumberMesh.hxx

compute_dof_pml, drop_interface_pml_dof

If the attribute compute_dof_pml is true, degrees of freedom associated with PMLs are numbered when NumberMesh is called. GetDofPML(i) will return the pml dof number of the i-th dof. If the i-th dof is not inside the PML, GetDofPML(i) returns -1. By default this attribute is set to false. If the attribute drop_interface_pml_dof is true, the dofs that are shared by the physical domain and PMLs will not be considered as PML. This approach will reduce the number of PML dofs.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// we consider a 2-D mesh numbering (with continuous elements)
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.compute_dof_pml = true; // you want to number dofs of PMLs
mesh_num.drop_interface_pml_dof = true; // dofs on the interface are not PML
mesh_num.NumberMesh();

// number of dofs in PMLs ?
int nb_dof_pml = mesh_num.GetNbDofPML();

Location :

NumberMesh.hxx

Init for Edge<Dimension2>

Syntax

void Init(int n1, int n2)
void Init(int n1, int n2, int ref)

You usually initialize the edge by specifying extremity numbers and optionally the reference number of the edge. This method only modifies the extremity numbers (and optionally the reference number), the other fields of the edge are not affected. We advise you to always sort numbers of the edge so that n1 is always lower than n2.

Example :

Mesh<Dimension2> mesh;

// first, vertices
mesh.ReallocateVertices(3);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0); mesh.Vertex(2) = R2(1,1);

// you specify two referenced edges
mesh.ReallocateBoundariesRef(2);

// then you can call Init
// the default reference number is 0
mesh.BoundaryRef(0).Init(0, 1);

// you can also specify a reference
mesh.BoundaryRef(1).Init(1, 2, 1);

// if you recall Init, it will keep previous reference (1 here)
mesh.BoundaryRef(1).Init(0, 2);

Location :

MeshElement.hxx
MeshElement.cxx

GetReference for Edge<Dimension2>

Syntax

int GetReference() const

This method returns the reference number associated with the edge. This reference number is usually determined in the definition of the geometry of the mesh. If you are using Gmsh, it is the number you specify just after Physical Line.

Example :

Mesh<Dimension2> mesh;
// you read a mesh
mesh.Read("test.mesh");
// then you loop over referenced edges
for (int i = 0: i < mesh.GetNbBoundaryRef(); i++)
  if (mesh.BoundaryRef(i).GetReference() == 2)
    cout << "Edge " << mesh.BoundaryRef(i) << endl;

Location :

MeshElement.hxx
MeshElement.cxx

SetReference for Edge<Dimension2>

Syntax

void SetReference(int ref)

This method changes the reference number associated with the edge.

Example :

Mesh<Dimension2> mesh;
// you read a mesh
mesh.Read("test.mesh");
// then you loop over referenced edges to change ref 2 to ref 3
for (int i = 0: i < mesh.GetNbBoundaryRef(); i++)
  if (mesh.BoundaryRef(i).GetReference() == 2)
    mesh.BoundaryRef(i).SetReference(3);

Location :

MeshElement.hxx
MeshElement.cxx

numVertex for Edge<Dimension2>

Syntax

int numVertex(int j) const

This method returns the vertex number of an extremity (first or second). You cannot change the vertex number, only retrieve it with that method. If you want to change extremities, you have to use Init.

Example :

Mesh<Dimension2> mesh;
// you read a mesh
mesh.Read("test.mesh");
// then you display extremities of referenced edges
for (int i = 0: i < mesh.GetNbBoundaryRef(); i++)
  {
    cout << "First extremity " << mesh.BoundaryRef(i).numVertex(0) << endl;
    cout << "Second extremity " << mesh.BoundaryRef(i).numVertex(1) << endl;
  }

Location :

MeshElement.hxx
MeshElement.cxx

GetNbVertices() for Edge<Dimension2>

Syntax

int GetNbVertices() const

This method returns 2 since there are two vertices on an edge.

Location :

MeshElement.hxx
MeshElement.cxx

ClearConnectivity for Edge<Dimension2>

Syntax

void ClearConnectivity()

This method removes informations about elements sharing the edge. After that call, there is no element associated with the edge. Usually, this low-level method is not used, since you can clear connectivity of all the mesh with method ClearConnectivity.

Example :

Mesh<Dimension2> mesh;
// you read a mesh
mesh.Read("test.mesh");

// we clear the elements shared by an edge ...
mesh.BoundaryRef(3).ClearConnectivity();

// if you want to clear the connectivity of all the mesh :
mesh.ClearConnectivity();

Location :

MeshElement.hxx
MeshElement.cxx

GetNbElements for Edge<Dimension2>

Syntax

int GetNbElements() const

This method returns the number of elements sharing the edge. If the edge is on the boundary, it will return 1 (only one element), otherwise it should return 2.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");
// special treatment for inside edges ?
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  if (mesh.Boundary(i).GetNbElements() == 2)
    cout << "Inside edge " << i << endl;

Location :

MeshElement.hxx
MeshElement.cxx

numElement for Edge<Dimension2>

Syntax

int numElement(int num) const

This methods returns the element number of an element adjacent to the edge.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");
// loop over referenced edges
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  {
    // element number
    int num_elem = mesh.BoundaryRef(i).numElement(0);
    // you can get global number of the edge
    int num_edge = i;
    // and local position of the edge on the element
    int num_loc = mesh.Element(num_elem).GetPositionBoundary(num_edge);
    cout << num_loc << "-th edge of element " << num_elem << endl;
  }

Location :

MeshElement.hxx
MeshElement.cxx

AddElement for Edge<Dimension2>

Syntax

void AddElement(int nf)

This methods adds an element to the edge. The number of elements is lower or equal to 2. If you try to add three elements to an edge, third element number will overwrite the second element number. This low level method should be used with caution. The research of adjacent elements for all the edges of the mesh is performed in method FindConnectivity (which is automatically called when reading/creating a mesh).

Example :

Mesh<Dimension2> mesh;

// first, vertices
mesh.ReallocateVertices(3);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0); mesh.Vertex(2) = R2(1,1);

// you specify two referenced edges
mesh.ReallocateBoundariesRef(2);
mesh.BoundaryRef(0).Init(0, 1, 1);
mesh.BoundaryRef(1).Init(1, 2, 1);

// then you have elements
mesh.ReallocateElements(1);
mesh.Element(0).InitTriangular(0, 1, 2, 1);

// you can specify elements by hand (not advised)
mesh.BoundaryRef(0).AddElement(0);
mesh.BoundaryRef(1).AddElement(0);

Location :

MeshElement.hxx
MeshElement.cxx

SetElement for Edge<Dimension2>

Syntax

void SetElement(int n0, int n1)

This methods sets the element numbers of an element adjacent to the edge. This low level method should be used with caution. The research of adjacent elements for all the edges of the mesh is performed in method FindConnectivity (which is automatically called when reading/creating a mesh).

Example :

Mesh<Dimension2> mesh;

// first, vertices
mesh.ReallocateVertices(3);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0); mesh.Vertex(2) = R2(1,1);

// you specify two referenced edges
mesh.ReallocateBoundariesRef(2);
mesh.BoundaryRef(0).Init(0, 1, 1);
mesh.BoundaryRef(1).Init(1, 2, 1);

// then you have elements
mesh.ReallocateElements(1);
mesh.Element(0).InitTriangular(0, 1, 2, 1);

// you can specify elements by hand (not advised)
mesh.BoundaryRef(0).AddElement(0);
mesh.BoundaryRef(1).AddElement(0);

// you can also provide the two elements 
mesh.BoundaryRef(0).SetElement(0, 3); // (element 0 and 3 are adjacent to this edge)

Location :

MeshElement.hxx
MeshElement.cxx

GetGeometryReference for Edge<Dimension2>

Syntax

int GetGeometryReference() const

This method returns the geometry number associated with the edge. This geometry number is only defined in .msh files generated by Gmsh, and corresponds to the identification number given to each line, circle and other curve in the .geo file. If you are reading meshes from other formats, the geometric reference will be 0 for every edge.

Example :

Mesh<Dimension2> mesh;
// you read a mesh generated by Gmsh
mesh.Read("test.msh");

// then you loop over referenced edges
for (int i = 0: i < mesh.GetNbBoundaryRef(); i++)
  if (mesh.BoundaryRef(i).GetGeometryReference() == 2)
    cout << "Edge " << mesh.BoundaryRef(i) << endl;

Location :

MeshElement.hxx
MeshElement.cxx

SetGeometryReference for Edge<Dimension2>

Syntax

void SetGeometryReference(int ref)

This method changes the geometry number associated with the edge.

Example :

Mesh<Dimension2> mesh;
// you read a mesh
mesh.Read("test.msh");
// then you loop over referenced edges to change geometry ref 2 to geometry ref 3
for (int i = 0: i < mesh.GetNbBoundaryRef(); i++)
  if (mesh.BoundaryRef(i).GetGeometryReference() == 2)
    mesh.BoundaryRef(i).SetGeometryReference(3);

Location :

MeshElement.hxx
MeshElement.cxx

IsTriangular for Face<Dimension2>

Syntax

bool IsTriangular() const

This method returns true if the face is a triangle, false otherwise.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
// specific treatment for triangles
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsTriangular())
    cout << "triangular element " << endl;

Location :

MeshElement.hxx
MeshElement.cxx

IsQuadrangular for Face<Dimension2>

Syntax

bool IsQuadrangular() const

This method returns true if the face is a quadrangle, false otherwise.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
// specific treatment for quadrangles
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsQuadrangular())
    cout << "quadrangular element " << endl;

Location :

MeshElement.hxx
MeshElement.cxx

GetHybridType for Face<Dimension2>

Syntax

int GetHybridType() const

This method returns 0 for a triangle, 1 for a quadrangle.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
// specific treatment for triangles/quadrangles
for (int i = 0; i < mesh.GetNbElt(); i++)
  switch (mesh.Element(i).GetHybridType())
    {
    case 0 :
      cout << "triangular element " << endl;
      break;
    case 1 :
      cout << "quadrangular element " << endl;
      break;
    }

Location :

MeshElement.hxx
MeshElement.cxx

GetOrientationEdge, GetOrientationBoundary for Face<Dimension2>

Syntax

bool GetOrientationEdge(int num_loc) const
int GetOrientationBoundary(int num_loc) const

This method returns the difference of orientation between local and global edge. It returns true if the orientations are equal, false if the orientations are opposite.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
int num_elem = 2;
for (int num_loc = 0; num_loc < mesh.Element(num_elem).GetNbEdges(); num_loc++)
  cout << "Difference of orientation " << mesh.Element(num_elem).GetOrientationEdge(num_loc);

Location :

MeshElement.hxx
MeshElement.cxx

numVertex for Face<Dimension2>

Syntax

int numVertex(int num_loc) const

This method returns the vertex number of a vertex of the face.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    int nb_vertices = mesh.Element(i).GetNbVertices();
    // storing vertices of the face in an array
    IVect num(nb_vertices);
    for (int j = 0; j < nb_vertices; j++)
      num(j) = mesh.Element(i).numVertex(j);
  }  

Location :

MeshElement.hxx
MeshElement.cxx

numEdge, numBoundary for Face<Dimension2>

Syntax

int numEdge(int num_loc) const
int numBoundary(int num_loc) const

This method returns the edge number of an edge of the face.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    int nb_edges = mesh.Element(i).GetNbEdges();
    // storing edges of the face in an array
    IVect num(nb_edges);
    for (int j = 0; j < nb_edges; j++)
      num(j) = mesh.Element(i).numEdge(j);
  }  

Location :

MeshElement.hxx
MeshElement.cxx

GetReference for Face<Dimension2>

Syntax

int GetReference() const

This method returns the reference number of the face. Usually this reference number is determined during the construction of the mesh. For example, in Gmsh it is the number just after Physical Surface.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetReference() == 2)
    cout << "Face with ref 2 " << endl;

Location :

MeshElement.hxx
MeshElement.cxx

GetNbVertices for Face<Dimension2>

Syntax

int GetNbVertices() const

This method returns the number of vertices of the face, e.g. three for a triangle.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetNbVertices() == 3)
    cout << "Triangular face " << endl;

Location :

MeshElement.hxx
MeshElement.cxx

GetNbEdges, GetNbBoundary for Face<Dimension2>

Syntax

int GetNbEdges() const
int GetNbBoundary() const

This method returns the number of edges of the face, e.g. 3 for a triangle.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetNbEdges() == 3)
    cout << "Triangular face " << endl;

Location :

MeshElement.hxx
MeshElement.cxx

SetVertex for Face<Dimension2>

Syntax

void SetVertex(int num_loc, int num_glob)

This method sets the global vertex number of a local vertex of the face. Usually, it is preferable to specify all the vertices of the face in the same time by calling methods InitTriangular, InitQuadrangular or Init.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// if you want to change a global number
int old_num = 23; int new_num = 57;
for (int i = 0; i < mesh.GetNbElt(); i++)
  for (int j = 0; j < mesh.Element(i).GetNbVertices(); j++)
    if (mesh.Element(i).numVertex(j) == old_num)
      mesh.Element(i).SetVertex(j, new_num);

Location :

MeshElement.hxx
MeshElement.cxx

SetReference for Face<Dimension2>

Syntax

void SetReference(int ref)

This method sets the reference number of the face.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
// if you want to change a reference
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetReference() == 2)
    mesh.Element(i).SetReference(3);

Location :

MeshElement.hxx
MeshElement.cxx

InitTriangular for Face<Dimension2>

Syntax

void InitTriangular(int n0, int n1, int n2, int ref)

This method specifies the three vertices of the face (which becomes a triangle) and the associated reference number.

Example :

Mesh<Dimension2> mesh;

// construction of a basic mesh
mesh.ReallocateVertices(3);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0); mesh.Vertex(2) = R2(1,1);

// creating elements
mesh.ReallocateElements(1);
mesh.Element(0).InitTriangular(0, 1, 2, 32);

Location :

MeshElement.hxx
MeshElement.cxx

InitQuadrangular for Face<Dimension2>

Syntax

void InitQuadrangular(int n0, int n1, int n2, int n3, int ref)

This method specifies the four vertices of the face (which becomes a quadrangle) and the associated reference number.

Example :

Mesh<Dimension2> mesh;
// construction of a basic mesh
mesh.ReallocateVertices(4);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0);
mesh.Vertex(2) = R2(1,1);  mesh.Vertex(3) = R2(0,1);

// creating elements
mesh.ReallocateElements(1);
mesh.Element(0).InitQuadrangular(0, 1, 2, 3, 32);

Location :

MeshElement.hxx
MeshElement.cxx

Init for Face<Dimension2>

Syntax

void Init(const IVect& num, int ref)

This method specifies the vertices of the face and the reference number. If the length of the array num is equal to 3, the face will be a triangle, if it is equal to 4, it will be a quadrangle.

Example :

Mesh<Dimension2> mesh;

// construction of a basic mesh
mesh.ReallocateVertices(4);
mesh.Vertex(0) = R2(0,0); mesh.Vertex(1) = R2(1,0);
mesh.Vertex(2) = R2(1,1);  mesh.Vertex(3) = R2(0,1);

// we create elements
mesh.ReallocateElements(1);
IVect num(4); num.Fill();
mesh.Element(0).Init(num, 32);

Location :

MeshElement.hxx
MeshElement.cxx

SetEdge for Face<Dimension2>

Syntax

void SetEdge(int num_loc, int num_glob)

This method sets the global number of an edge of the face. This method is a low-level, since global edge numbers are attributed when calling FindConnectivity. In regular use, this method should not be called.

Example :

Mesh<Dimension2> mesh;
mesh.ReadMesh("test.mesh");

// changing first edge number of element 3 (not advised to do that)
mesh.Element(3).SetEdge(0, 32);

Location :

MeshElement.hxx
MeshElement.cxx

GetPositionBoundary for Face<Dimension2>

Syntax

int GetPositionBoundary(int num_glob) const

This method returns the local position of a global edge.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// loop over referenced edges
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  {
    // element number
    int num_elem = mesh.BoundaryRef(i).numElement(0);
    // usually referenced edges are placed before the other edges
    int num_edge = i;
    // local position of the edge on the element
    int num_loc = mesh.Element(num_elem).GetPositionBoundary(num_edge);
    cout << num_loc << "-th edge of element " << num_elem << endl;
}

Location :

MeshElement.hxx
MeshElement.cxx

ClearConnectivity for Face<Dimension2>

Syntax

void ClearConnectivity()

This method clears the edges which have been found. The edge numbers are set to -1 after a call to that function. This low-level method should not be used, you can clear connectivity of all the mesh by calling ClearConnectivity.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");

// and after you can clear connectivities of elements only
for (int i = 0; i < mesh.GetNbElt(); i++)
  mesh.Element(i).ClearConnectivity();

Location :

MeshElement.hxx
MeshElement.cxx

ClearEdges for Face<Dimension2>

Syntax

void ClearEdges()

This method clears the edges which have been found. The edge numbers are set to -1 after a call to that function. This low-level method should not be used, you can clear connectivity of all the mesh by calling ClearConnectivity.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// and after you clear informations about edges
for (int i = 0; i < mesh.GetNbElt(); i++)
  mesh.Element(i).ClearEdges();

Location :

MeshElement.hxx
MeshElement.cxx

IsPML for Face<Dimension2>

Syntax

bool IsPML() const

This method returns true if the face is considered as a face in PML layers. This status can be modified with methods SetPML() and UnsetPML().

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// then you can count pml elements for example
int nb_elt_pml = 0;
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsPML())
    nb_elt_pml++;

Location :

MeshElement.hxx
MeshElement.cxx

SetPML for Face<Dimension2>

Syntax

void SetPML(int num, int type)

This method informs that the element is in a PML. The first argument num if the PML region number, whereas the second argument type is the type of PML.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// you can mark some elements as pml
mesh.Element(23).SetPML(0, 0);

Location :

MeshElement.hxx
MeshElement.cxx

UnsetPML for Face<Dimension2>

Syntax

void UnsetPML()

This method informs that the element is not in a PML.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");
// you can mark some elements as pml
mesh.Element(23).SetPML(0, 0);
// and unmark them
mesh.Element(23).UnsetPML();

Location :

MeshElement.hxx
MeshElement.cxx

GetMemorySize for Face<Dimension2>

Syntax

size_t GetMemorySize() const

This method returns the size used by a face in bytes.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// loop over elements
for (int i = 0; i < mesh.GetNbElt(); i++)
  cout << mesh.Element(i).GetMemorySize() << endl;

Location :

MeshElement.hxx
MeshElement.cxx

GetTypePML for Face<Dimension2>

Syntax

int GetTypePML() const

This method returns the type of PML associated with the face.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// loop over PML elements
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsPML())
    cout << mesh.Element(i).GetTypePML() << endl;

Location :

MeshElement.hxx
MeshElement.cxx

GetNumberPML for Face<Dimension2>

Syntax

int GetNumberPML() const

This method returns the PML region number associated with the face.

Example :

Mesh<Dimension2> mesh;

// we read a mesh
mesh.Read("test.mesh");

// loop over PML elements
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsPML())
    cout << mesh.Element(i).GetNumberPML() << endl;

Location :

MeshElement.hxx
MeshElement.cxx

IsCurved for Face<Dimension2>

Syntax

bool IsCurved() const

This method returns true if the face is considered as a curved face. This status can be modified with methods SetCurved() and UnsetCurved().

Example :

Mesh<Dimension2> mesh;
// we read a mesh 
// if the mesh file is appropriate, there may be curved elements
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// then you can count curved elements
int nb_elt_curved = 0;
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).IsCurved())
    nb_elt_curved++;

Location :

MeshElement.hxx
MeshElement.cxx

SetCurved for Face<Dimension2>

Syntax

void SetCurved()

This method informs that the element is curved.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");
// you can mark some elements as curved
mesh.Element(23).SetCurved();

Location :

MeshElement.hxx
MeshElement.cxx

UnsetCurved for Face<Dimension2>

Syntax

void UnsetCurved()

This method informs that the element is not curved.

Example :

Mesh<Dimension2> mesh;
// we read a mesh
mesh.Read("test.mesh");
// you can mark some elements as curved
mesh.Element(23).SetCurved();
// and unmark them
mesh.Element(23).UnsetCurved();

Location :

MeshElement.hxx
MeshElement.cxx

Boundary

Syntax

Edge<Dimension2>& Boundary(int i)

This method gives access to an edge of the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we can browse all the edges of the mesh,
// to find an edge of extremities 0 and 2 for example
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  {
    int n1 = mesh.Boundary(i).numVertex(0);
    int n2 = mesh.Boundary(i).numVertex(1);
    if ((n1 == 0)&&(n2 == 2))
      cout << "We found the edge 0-2" << endl;
  }

Related topics :

GetNbBoundary
FindConnectivity

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetNbBoundary

Syntax

int GetNbBoundary() const

This method returns the number of edges present in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we print the number of edges
cout << "The mesh contains " << mesh.GetNbBoundary() << "edges" << endl;

Related topics :

Boundary
FindConnectivity

Location :

Mesh2D.hxx
Mesh2DInline.cxx

ReallocateBoundariesRef

Syntax

void ReallocateBoundariesRef(int n)

This method reallocates the array containing the referenced edges. Notice that previous edges may be lost during the process. This method is usually called when a new mesh is created (for example a mesh of a rectangle).

Example :

Mesh<Dimension2> mesh;
// we assume that vertices are constructed
// you can specify two referenced edges
mesh.ReallocateBoundariesRef(2);
int ref = 5;
mesh.BoundaryRef(0).Init(0, 1, ref);
mesh.BoundaryRef(1).Init(2, 3, ref);

Related topics :

ResizeBoundariesRef

Location :

Mesh2D.hxx
Mesh2DInline.cxx

ResizeBoundariesRef

Syntax

void ResizeBoundariesRef(int n)

This method reallocates the array containing the referenced edges. Notice that previous edges are conserved during the process (edges not beyond the upper bound of the array). This method is usually called when an operation on the mesh is performed (like splitting the mesh).

Example :

Mesh<Dimension2> mesh;
// we assume that vertices are constructed
// you can specify two referenced edges
mesh.ReallocateBoundariesRef(2);
int ref = 5;
mesh.BoundaryRef(0).Init(0, 1, ref);
mesh.BoundaryRef(1).Init(2, 3, ref);
// and change the size while keeping the two previous entries
mesh.ResizeBoundariesRef(3);
// and you specify only the new referenced edges
mesh.BoundaryRef(2).Init(1, 3, ref);

Related topics :

ReallocateBoundariesRef

Location :

Mesh2D.hxx
Mesh2DInline.cxx

BoundaryRef

Syntax

Edge<Dimension2>& BoundaryRef(int i)

This method gives access to a referenced edge of the mesh (usually edges on a specific boundary of the mesh). These edges are present after the reading of a mesh. If you modify a referenced edge manually, you should sort the numbers, so that the first extremity number is lower than the second one. If you have provided unsorted numbers, please call the method SortBoundariesRef. A modification of referenced edges should be followed by a call to ConstructCrossReferenceBoundaryRef or FindConnectivity. If you have curved referenced edges, you may have to call SetPointInsideEdge in order to specify position of internal nodes for each referenced edge.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we can browse all the referenced edges of the mesh,
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  {
    int ref = mesh.BoundaryRef(i).GetReference();
    cout << "Reference of the edge : " << ref << endl;
  }

Related topics :

GetNbBoundaryRef

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetNbBoundaryRef

Syntax

int GetNbBoundaryRef() const
int GetNbEdgesRef() const

This method returns the number of referenced edges present in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we print the number of referenced edges
cout << "The mesh contains " << mesh.GetNbBoundaryRef() << "referenced edges" << endl;

Related topics :

BoundaryRef

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetTriangleReferenceElement

Syntax

const TriangleGeomReference& GetTriangleReferenceElement() const

This method returns the nodal triangle used for curved triangles.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.msh");
// and we retrieve the nodal triangle
const TriangleGeomReference& tri = mesh.GetTriangleReferenceElement();
cout << "Number of nodal points on the triangle " << tri.GetNbPointsNodalElt() << endl;

Related topics :

TriangleGeomReference

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetQuadrilateralReferenceElement

Syntax

const QuadrangleGeomReference& GetQuadrilateralReferenceElement() const

This method returns the nodal quadrilateral used for curved quadrilaterals.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(2);
mesh.Read("test.msh");
// and we retrieve the nodal quadrangle
const QuadrangleGeomReference& quad = mesh.GetQuadrilateralReferenceElement();
cout << "Number of nodal points on the quadrangle " << quad.GetNbPointsNodalElt() << endl;

Related topics :

QuadrangleGeomReference

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetReferenceElementSurface

Syntax

void GetReferenceElementSurface(Vector<const ElementGeomReference<Dimension1>* >& elt) const

This method fills the array elt with the 1-D finite element used for curved boundaries. The length of the array elt is equal to 1.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.msh");

// and we retrieve the 1-D nodal finite element
Vector<const ElementGeomReference<Dimension1>* > elt1D;
mesh.GetReferenceElementSurface(elt1D);
cout << "Nodal 1D points " << elt1D(0)->Points() << endl;

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetReferenceElementVolume

Syntax

void GetReferenceElementVolume(Vector<const ElementGeomReference<Dimension2>* >& elt) const

This method fills the array elt with the 2-D finite element used for each element.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.msh");

// and we retrieve the 2-D nodal finite elements
Vector<const ElementGeomReference<Dimension2>* > elt2D;
mesh.GetReferenceElementVolume(elt2D);
cout << "number of nodal points on element 5 = " << elt2D(5)->GetNbPointsNodalElt() << endl;

Location :

Mesh2D.hxx
Mesh2D.cxx

GetLocalNodalNumber

Syntax

Vector<int> GetLocalNodalNumber(int type, int num_loc)

This method returns the nodal numbers of the edge num_loc of a triangle (if type is 0) or a quadrangle (if type is 1).

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.msh");

// number of nodal points for the edge 2 of the unit quadrangle
int num_loc = 2;
IVect num_nodes = mesh.GetLocalNodalNumber(1, num_loc);

Location :

Mesh2D.hxx
Mesh2D.cxx

HybridMesh

Syntax

bool HybridMesh() const

This method returns true if the mesh contains both quadrilaterals and triangles. If the mesh contains only triangles or only quadrilaterals, it will return false.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// it is an hybrid mesh ?
if (mesh.HybridMesh())
  cout << "The mesh is hybrid " << endl;

Related topics :

IsOnlyTriangular
IsOnlyQuadrilateral

Location :

Mesh2D.hxx
Mesh2D.cxx

IsOnlyTriangular

Syntax

bool IsOnlyTriangular() const

This method returns true if the mesh contains only triangles. If the mesh contains some quadrilaterals, it will return false.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// is it a pure triangular mesh ?
if (mesh.IsOnlyTriangular())
  cout << "The mesh contains only triangles " << endl;

Related topics :

HybridMesh
IsOnlyQuadrilateral

Location :

Mesh2D.hxx
Mesh2DInline.cxx

IsOnlyQuadrilateral

Syntax

bool IsOnlyQuadrilateral() const

This method returns true if the mesh contains only quadrilaterals. If the mesh contains some triangles, it will return false.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// is it a pure quadrilateral mesh ?
if (mesh.IsOnlyQuadrilateral())
  cout << "The mesh contains only quadrilaterals " << endl;

Related topics :

HybridMesh
IsOnlyTriangular

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetMemorySize

Syntax

size_t GetMemorySize() const

This method returns the memory used to store the mesh in bytes.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.msh");

// you want to know the memory used to store the object mesh :
cout << "Memory used to store mesh = " << mesh.GetMemorySize() << " bytes " << endl;

Location :

Mesh2D.hxx
Mesh2D.cxx

GetNbTriangles

Syntax

int GetNbTriangles() const

This method returns the number of triangles present in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we print the number of triangles
cout << "The mesh contains " << mesh.GetNbTriangles() << " triangles" << endl;

Related topics :

GetNbQuadrangles

Location :

Mesh2D.hxx
Mesh2DInline.cxx

GetNbQuadrangles

Syntax

int GetNbQuadrangles() const

This method returns the number of quadrilaterals present in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// and we print the number of quadrangles
cout << "The mesh contains " << mesh.GetNbQuadrangles() << " quadrangles" << endl;

Related topics :

GetNbTriangles

Location :

Mesh2D.hxx
Mesh2DInline.cxx

SetGeometryOrder

Syntax

void SetGeometryOrder(int r)
void SetGeometryOrder(int r, bool only_quad)

This method changes the order of approximation used for the definition of the geometry. It is necessary to call that method before reading the mesh, otherwise the mesh will be converted into a first order mesh, and you will lose the information about curved boundaries contained in the mesh file. Usually in Montjoie codes (except those in the directory src/Program/Mesh), this step is performed when the data file is read (for the field OrderDiscretization). When a mesh is read, it is translated into Montjoie structures, so at the order specified when you have called SetGeometryOrder (the default order is 1). You can change the order of approximation after reading the mesh, in that case it will perform an interpolation to obtain the new nodal points on referenced edges of the mesh or exploit the informations provided for the curves. Of course, if you have decreased the order of approximation, some information is lost, and if you get back to the initial order of approximation, the internal points will not be exactly the same. The second argument is optional, if you put true, only quadrilateral finite element is computed. It is assumed that the mesh will not contain triangles in that case.

Example :

Mesh<Dimension2> mesh;
// first you specify which order you wish
mesh.SetGeometryOrder(3);
// you can specify a given curve if you want (you don't need if the input mesh is already curved)
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
// if you don't provide parameters with method SetCurveParameters, the
// parameters will be automatically found 
// we read the mesh
mesh.Read("test.mesh");
// you can decrease the order of approximation if you want
mesh.SetGeometryOrder(2);

Related topics :

Read

Location :

Mesh2D.hxx
Mesh2D.cxx

GetUpperBoundNumberOfEdges

Syntax

int GetUpperBoundNumberOfEdges() const

This method returns an upper bound of the number of edges. It is used by the method FindConnectivity, otherwise the user should not call that method.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// upper bound of the number of edges
int nb_max_edges = mesh.GetUpperBoundNumberOfEdges();
// then you can compare with the exact number
int nb_edges = mesh.GetNbEdges();

Related topics :

GetNbEdges

Location :

Mesh2D.hxx
Mesh2D.cxx

GetNodesCurvedMesh

Syntax

int GetNodesCurvedMesh(VectR2& PosNodes, Vector<IVect>& Nodle, int rmax) const
int GetNodesCurvedMesh(VectR2& PosNodes, Vector<IVect>& Nodle, int rmax, bool lobatto) const

This method computes the positions of each node of the mesh and the array containing node numbers for each element. For first-order meshes, the nodes of the mesh are the vertices, and Nodle will contain the vertex numbers for each element. For high-order meshes, vertices and intern nodes are computed, and nodes numbers are stored for each element. The local numbering of nodes on a triangle is displayed in the method ConstructTriangularNumbering, and on a quadrilateral, you can look at the picture of the method ConstructQuadrilateralNumbering. You can require to have nodes for an order rmax lower than the order of approximation of the mesh, but if you ask nodes for order rmax greater, it will give you nodes for the order of approximation of the mesh, and the method will return this order. If you want to avoid this drawback, you can call SetGeometryOrder before. If the last argument is set to true, the nodal points are computed for Gauss-Lobatto points. Otherwise, the nodal points are computed for a regular subdivision.

Example :

Mesh<Dimension2> mesh;
mesh.SetGeometryOrder(3);
// we read the mesh
mesh.Read("test.mesh");
// then you want to explicit all the nodes and node numbering
VectR2 PosNodes; Vector<IVect> Nodle;
// here you will get only third-order mesh (r = 3)
int r = mesh.GetNodesCurvedMesh(PosNodes, Nodle, 4);
cout << "Actual order of nodes " << r << endl;
// printing all the nodes
cout << "Position of nodes " << endl;
for (int i = 0; i < PosNodes.GetM(); i++)
  cout << PosNodes(i) << endl;
    
for (int i = 0; i < Nodle.GetM(); i++)
  cout << "Nodes of element " << i << endl << Nodle(i) << endl;
   

Related topics :

SetGeometryOrder

Location :

Mesh2D.hxx
Mesh2D.cxx

Read

Syntax

void Read(string file_name)

This method reads a mesh from a file. The mesh format is automatically found with the extension of the file name. If the mesh file is curved, you will have to call SetGeometryOrder before, so that the mesh is converted to an high-order mesh with the order specified by the call of SetGeometryOrder. In Montjoie, we use generally isoparametric elements, the order of approximation is usually the same for the geometry and the finite element space.

Example :

Mesh<Dimension2> mesh;
mesh.SetGeometryOrder(3);
// we read the mesh
mesh.Read("test.mesh");
   

Related topics :

SetGeometryOrder

Location :

Mesh2D.hxx
Mesh2D.cxx

Write

Syntax

void Write(string file_name) const

This method writes a mesh in a file. The mesh format is automatically found with the extension of the file name.

Example :

Mesh<Dimension2> mesh;
mesh.SetGeometryOrder(3);
// we read the mesh
mesh.Read("test.mesh");
// and writing in another format
mesh.Write("test.msh");

Related topics :

Read

Location :

Mesh2D.hxx
Mesh2D.cxx

ReorientElements

Syntax

void ReorientElements()

This method changes the local numbering of elements, so that the jacobian is positive for all elements. This method is called when a mesh is read from a file so that you probably don't need to call it. For most of the methods (except low-level methods such as Vertex, Element, etc) defined in class Mesh, the orientation of elements is set so that jacobians are always positive.

Example :

Mesh<Dimension2> mesh;
// if you construct by your-self a mesh
mesh.ReallocateVertices(4);
mesh.Vertex(0).Init(-1.0, -1.0);
// etc
mesh.ReallocateElements(5);
mesh.Element(0).InitQuadrangular(4, 2, 5, 8, 1);
// etc

// you will need to reorient elements such that jacobian is positive
mesh.ReorientElements();

Related topics :

Read

Location :

Mesh2D.hxx
Mesh2D.cxx

FindConnectivity

Syntax

void FindConnectivity()

This method finds all the edges of the mesh (we may call them global edges). Referenced edges are constituting the first global edges (internal edges are placed at the end of the array containing all the edges).

Example :

Mesh<Dimension2> mesh;

// we read the mesh (here FindConnectivity is called by Read)
mesh.Read("test.mesh");

// you can ask to recompute edges (not needed in this case)
mesh.FindConnectivity()

// printing the number of edges
cout << "Number of edges " << mesh.GetNbEdges() << endl;

Related topics :

ClearConnectivity

Location :

Mesh2D.hxx
Mesh2D.cxx

RemoveDuplicateVertices

Syntax

void RemoveDuplicateVertices()

This method removes from the mesh duplicate vertices, i.e. vertices which are twice in the mesh, so that in the new mesh, each vertex is unique.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we remove duplicate vertices
mesh.RemoveDuplicateVertices();
// printing the new number of vertices
cout << "Number of vertices " << mesh.GetNbVertices() << endl;

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

ForceCoherenceMesh

Syntax

void ForceCoherenceMesh()
void ForceCoherenceMesh(bool keep_edge)

This method removes from the mesh unused vertices, i.e. vertices which are not belonging to elements. The referenced edges that do not belong to elements are removed. If the optional argument keep_edge is equal to true, the referenced edges are not removed, and vertices that belong to reference edges are also conserved.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we remove unused vertices
mesh.ForceCoherenceMesh();
// printing the new number of vertices
cout << "Number of vertices " << mesh.GetNbVertices() << endl;

mesh.Read("test.mesh");
// in that case, we keep edges and vertices of edges (not only vertices of elements)
mesh.ForceCoherenceMesh(true);
// printing the new number of vertices
cout << "Number of vertices " << mesh.GetNbVertices() << endl;

Location :

Mesh2D.hxx
Mesh2D.cxx

ClearZeroBoundaryRef

Syntax

void ClearZeroBoundaryRef()

This method removes from the mesh referenced edges whose reference is equal to 0.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we remove edges of reference 0
mesh.ClearZeroBoundaryRef();
// printing the new number of referenced edges
cout << "Number of reference edges " << mesh.GetNbBoundaryRef() << endl;

Location :

Mesh2D.hxx
Mesh2D.cxx

ApplyLocalRefinement

Syntax

void ApplyLocalRefinement(int n, int nb_levels, Real_wp ratio)

This method performs a local refinement around the n-th vertex. nb_levels is the number of refinement steps performed, ratio is the geometric progression used. A ratio equal to 2 is a standard value, since in that case each edge adjacent to the vertex is divided into two equal parts. For ratio different from 2, each edge is separated into unequal parts of relative length 1/ratio (near the vertex) and 1 - 1/ratio.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we apply a local refinement near vertex 32
mesh.ApplyLocalRefinement(32, 3, 2.0);

In the figures below, we show how a triangular and quadrilateral mesh is modified when calling this method with a ratio equal to 3.0

Loading     ⇒     Loading
Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

AddVertexConformal

Syntax

void AddVertexConformal(R2 point)

This method adds a new vertex in the mesh. This new vertex is linked with old vertices such that the mesh stays conformal. Usually, it implies the creation of triangles.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// we add a vertex in the mesh
mesh.AddVertexConformal(R2(4.0, -0.8));

Location :

Mesh2D.hxx
Mesh2D.cxx

SymmetrizeMesh

Syntax

void SymmetrizeMesh(int num_coor, Real_wp xc, bool keep_edge)

This method symmetrizes the mesh with respect to coordinate num_coor (0 : x, 1 : y). The axis of symmetry is determined by the equation x = xc (or y = xc if num_coor is equal to 1). By default, the referenced edges present on the symmetry axis are removed, but you can ask to keep them.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we symmetrize mesh with respect to axis y = -2
mesh.SymmetrizeMesh(1, -2.0);

In the figure, we show the effect of this method on a triangular mesh.

Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

PeriodizeMeshTeta

Syntax

void PeriodizeMeshTeta(const R2& center)

This method generates a complete cyclic mesh with an initial cell. You specify the center of the polar coordinates, usually (0, 0).

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we symmetrize mesh by cyclicity
R2 center;
mesh.PeriodizeMeshTeta(center);

In the figure, we show the effect of this method on a triangular mesh.

Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

SplitIntoQuadrilaterals

Syntax

void SplitIntoQuadrilaterals()

This method transforms a triangular or hybrid mesh into a purely quadriteral mesh by splitting each triangle into three quadrilaterals, and each quadrangle into four quadrilaterals as shown in the figure below.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we split initial mesh into quadrilaterals only
mesh.SplitIntoQuadrilaterals();
Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

SplitIntoTriangles

Syntax

void SplitIntoTriangles()

This method transforms a quadrilateral or hybrid mesh into a purely triangular mesh by splitting each quadrilateral into two triangles as shown in the figure below.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we split initial mesh into triangles only
mesh.SplitIntoTriangles();
Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

GetMeshSubdivision

Syntax

void GetMeshSubdivision(const VectReal_wp& outside_subdiv, const Vector<Vector<R2> >& points_div, Vector<R2>& PosNodes, Vector<IVect>& Nodle) const

This method returns a subdivision of the mesh through arrays PosNodes (nodal points) and Nodle (connectivity). Local numbering of triangle and quadrilateral is explained in method ConstructTriangularNumbering and ConstructQuadrilateralNumbering.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// defining regular subdivision
int order = 3;
VectReal_wp step_x(order+1);
for (int i = 0; i ≤ order; i++)
  step_x(i) = Real_wp(i)/order;
	
Vector<VectR2> points_div(2);
points_div(0).Reallocate((order+1)*(order+2)/2);
points_div(1).Reallocate((order+1)*(order+1));

// regular points on triangles
Matrix<int> NumNodes2D, coor;
MeshNumbering<Dimension2>::ConstructTriangularNumbering(order, NumNodes2D, coor);
for (int i = 0; i ≤ order; i++)
   for (int j = 0; j ≤ order-i; j++)
      points_div(0)(NumNodes2D(i, j)).Init(step_x(i), step_x(j));
	
// regular points on quadrangles
MeshNumbering<Dimension2>::ConstructQuadrilateralNumbering(order, NumNodes2D, coor);
for (int i = 0; i ≤ order; i++)
   for (int j = 0; j ≤ order; j++)
      points_div(1)(NumNodes2D(i, j)).Init(step_x(i), step_x(j));
	
// we get a subdivision of the mesh
VectR2 PosNodes; Vector<IVect> Nodle;
mesh.GetMeshSubdivision(step_x, points_div, PosNodes, Nodle);

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

SubdivideMesh

Syntax

void SubdivideMesh(const VectReal_wp& outside_subdiv)
void SubdivideMesh(const VectReal_wp& outside_subdiv, const Vector<VectR2>& points_surf, Vector<VectR2>& points_div, Vector<IVect>& Nodle, Vector<IVect>& NodleSurf)

This method subdivides a mesh by giving subdivision on the edge [0,1] in array outside_subdiv. You can retrieve the subdivided points of triangle and quadrangle in array points_div, and the new vertex numbers on each element in array Nodle. Local numbering of triangle and quadrilateral is explained in method ConstructTriangularNumbering and ConstructQuadrilateralNumbering.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// defining regular subdivision
int order = 3;
VectReal_wp step_x(order+1);
for (int i = 0; i ≤ order; i++)
  step_x(i) = Real_wp(i)/order;

// subdividing mesh
mesh.SubdivideMesh(step_x);	

We see below the result obtained on a mesh by this method

Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

SymmetrizeMeshToGetPositiveX

Syntax

void SymmetrizeMeshToGetPositiveX()

This method performs a symmetry of the mesh with respect to axis x = 0, if the mesh is initially placed in area x < 0. This is used for axisymmetric computation to enforce a mesh placed in area x ≥ 0. The transformation of a mesh is shown in the figure below.

Loading     ⇒     Loading

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// transforming mesh
mesh.SymmetrizeMeshToGetPositiveX();

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

CreateRegularMesh

Syntax

void CreateRegularMesh(R2 ptMin, R2 ptMax, TinyVector<int,2> nbPoints, int ref_domain, TinyVector<int,4> ref_boundary, int type_mesh, R2 ratio)

This method creates a regular mesh of a rectangular zone defined by the left bottom corner ptMin and the right top corner ptMax. You specify the number of points of the mesh on each direction (x and y) in nbPoints, the reference to attribute to each boundary is given in ref_boundary. type_mesh accepts TRIANGULAR_MESH, QUADRILATERAL_MESH, RADAU_MESH, so that you can construct purely triangular or quadriteral regular mesh. The last argument ratio is optional. If present, the positions will follow a geometric progression (ratio(0) is the common ratio for x-coordinate, ratio(1) for y-coordinate).

Parameters

ptMin (in)
left bottom corner of the rectangle to mesh
ptMax (in)
right top corner of the rectangle to mesh
nbPoints (in)
number of points in x and y
ref_domain (in)
reference of created elements
ref_boundary (in)
reference of boundary edges
type_mesh (in)
type of mesh (TRIANGULAR_MESH or QUADRILATERAL_MESH)
ratio (optional)
geometric progression in x and y

Example :

Mesh<Dimension2> mesh;
// we construct a regular mesh of square [-2, 2]^2
// with 20 points along x, 11 points along y, quadrilateral mesh
TinyVector<int, 4> ref_boundary;
ref_boundary.Fill(3);
mesh.CreateRegularMesh(R2(-2, -2), R2(2, 2), 1, TinyVector<int, 2>(20, 11),
                       ref_boundary, mesh.QUADRILATERAL_MESH);

Related topics :

FileMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

ExtrudeCoordinate

Syntax

void ExtrudeCoordinate(int num_coor, int nb_layers, Real_wp pos, Real_wp delta)

This method adds layers on a cartesian boundary of a mesh, it is actually used in Montjoie to add PML layers around the computational domain. num_coor is equal to 0 or 1, to add layers in x-coordinate or y-coordinate. Layers are then added to fill an area between x = pos and x = pos + delta (or y if num_coor is equal to 1). These layers comprise only quadrilaterals.

Parameters

num_coor (in)
direction of extrusion (0 : x, 1 : y)
nb_layers (in)
number of layers to add
pos (in)
extrusion starts from x=pos (or y=pos)
delta (in)
thickness of created layers

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// we add layers between x = -2 and x = -3
int nb_layers = 3;
mesh.ExtrudeCoordinate(0, nb_layers, -2.0, -1.0);

In the figure below, we show the effect of this method on a triangular mesh.

Loading     ⇒     Loading

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

ConstructMesh

Syntax

void ConstructMesh(int type_mesh, VectString param)

This method constructs a mesh with parameters given in array param. These parameters are the same as parameters described in field FileMesh. If type_mesh is equal to QUADRILATERAL_MESH, the mesh will be split into quadrilaterals if it is containing triangles. If type_mesh is equal to TRIANGULAR_MESH, the mesh will be split into triangles if it is containing quadrilaterals. If type_mesh is equal to HYBRID_MESH, no splitting is performed.

Example :

Mesh<Dimension2> mesh;
VectString param(1);
param(0) = string("carre.mesh");
mesh.ConstructMesh(mesh.HYBRID_MESH, param);

Related topics :

FileMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

RearrangeElements

Syntax

void RearrangeElements()

This method changes the ordring of elements, so that the triangles are placed first in the list of elements. The method ConstructMesh calls already this method.

Example :

  Mesh<Dimension2> mesh;
  mesh.Read("test.mesh");

  // if you want to be sure that triangles are placed first
  mesh.RearrangeElements();

  int N = mesh.GetNbTriangles();
  int nb_elt = mesh.GetNbElt();
  // triangles will be between 0 and N-1, quadrangles between N and nb_elt-1

Location :

Mesh2D.hxx
Mesh2D.cxx

CreateSubmesh

Syntax

void CreateSubmesh(Mesh<Dimension2>& sub_mesh, int nb_vertices, int nb_elt, Vector<bool> vertex_on_subdomain, Vector<bool> element_on_subdomain)

This method extracts a part of a mesh to create a small sub-mesh.

Parameters

sub_mesh (out)
created sub-mesh
nb_vertices (in)
number of vertices to extract
nb_elt (in)
number of elements to extract
vertex_on_subdomain (in)
if vertex_on_subdomain(i) is true the vertex i must be extracted
element_on_subdomain (in)
if element_on_subdomain(i) is true the element i must be extracted

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

// for example, we are marking elements with reference equal to 3
int nb_vertices = 0, nb_elt = 0;
Vector<bool> vertex_on_subdomain(mesh.GetNbVertices()), element_on_subdomain(mesh.GetNbElt());
vertex_on_subdomain.Fill(false);
element_on_subdomain.Fill(false);

for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetReference() == 3)
    {
      element_on_subdomain(i) = true;
      nb_elt++;
      for (int j = 0; j < mesh.Element(i).GetNbVertices(); j++)
        {
          int nv = mesh.Element(i).numVertex(j);
          if (!vertex_on_subdomain(nv))
            {
              vertex_on_subdomain(nv) = true;
              nb_vertices++;
            }
        }
    }

// and extracting the sub-mesh
Mesh<Dimension2> sub_mesh;
mesh.CreateSubmesh(sub_mesh, nb_vertices, nb_elt, vertex_on_subdomain, element_on_subdomain);

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

AppendMesh

Syntax

void AppendMesh(Mesh<Dimension2> mesh, bool elimine = true, R2 u = 0)

This method appends a mesh to the current mesh. Duplicate vertices are removed, except if you are providing a second argument equal to false. The last optional argument is a translation vector (the mesh is appended with a translation).

Parameters

mesh (in)
mesh to append
elimine (optional)
if false, duplicate vertices are kept
u (optional)
translation vector

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// reading a small mesh
mesh2.Read("small.mesh");
// and appending to the first mesh
mesh.AppendMesh(mesh2);

// if you want to append again the same mesh, but with a translation
mesh.AppendMesh(mesh2, false, R2(5, 0));

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

FjElemNodal

Syntax

void FjElemNodal(const VectR2& s, SetPoints<Dimension2>& Pts, const Mesh<Dimension2>& mesh, int n) const

This method computes nodal points of an element.

Parameters

s (in)
vertices of the element
Pts (out)
computed nodal points
mesh (in)
current mesh
n (in)
element number

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

for (int j = 0; j < Pts.GetNbPointsNodal(); j++)
  cout << "Nodal point " << j << " : " << Pts.GetPointNodal(j) << endl;

Location :

MeshBase.hxx
MeshBase.cxx

DFjElemNodal

Syntax

void DFjElemNodal(const VectR2& s, const SetPoints<Dimension2>& Pts, SetMatrices<Dimension2>& Mat, const Mesh<Dimension2>& mesh, int n) const

This method computes jacobian matrices on nodal points of an element.

Parameters

s (in)
vertices of the element
Pts (in)
nodal points
Mat (out)
computed jacobian matrices
mesh (in)
current mesh
n (in)
element number

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

// then jacobian matrices on these nodal points
SetMatrices<Dimension2> Pts;
mesh.DFjElemNodal(s, Pts, Mat, mesh, n);

Location :

MeshBase.hxx
MeshBase.cxx

Fj

Syntax

void Fj(const VectR2& s, const SetPoints<Dimension2>& Pts, const R2& pt_loc, R2& pt_glob, const Mesh<Dimension2>& mesh, int n) const

This method computes the image by transformation Fi of a local point placed on the reference element.

Parameters

s (in)
vertices of the element
Pts (in)
nodal points
pt_loc (in)
local coordinates of point
pt_glob (out)
global coordinates of point
mesh (in)
current mesh
n (in)
element number

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

// then computation of pt_glob = Fi(pt_loc)
R2 pt_loc(0.3, 0.2);
R2 pt_glob;
mesh.Fj(s, Pts, pt_loc, pt_glob, mesh, n);

Location :

MeshBase.hxx
MeshBase.cxx

DFj

Syntax

void DFj(const VectR2& s, const SetPoints<Dimension2>& Pts, const R2& pt_loc, Matrix2_2& dfj, const Mesh<Dimension2>& mesh, int n) const

This method computes jacobian matrix DFi on a point of the reference element.

Parameters

s (in)
vertices of the element
Pts (in)
nodal points
pt_loc (in)
local coordinates of point
dfj (out)
jacobian matrix
mesh (in)
current mesh
n (in)
element number

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

// then computation of dfj = DFi(pt_loc)
R2 pt_loc(0.3, 0.4);
Matrix2_2 dfj;
mesh.DFj(s, Pts, pt_loc, dfj, mesh, n);

Location :

MeshBase.hxx
MeshBase.cxx

GetNbPointsNodalElt

Syntax

int GetNbPointsNodalElt(int i)

This method returns the number of nodal points associated with the element i.

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// number of nodal points of element 18 ?
int nb_pts = mesh.GetNbPointsNodalElt(18);

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetNbPointsNodalBoundary

Syntax

int GetNbPointsNodalBoundary(int num_elem, int num_loc)

This method returns the number of nodal points associated with an edge.

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// number of nodal points of an edge
int nb_pts = mesh.GetNbPointsNodalBoundary(18, 0);

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetNodalNumber

Syntax

int GetNodalNumber(int num_elem, int num_loc, int k)

This method returns the number of the k-th nodal point of the edge num_loc of element num_elem.

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// number of a nodal point ?
int num_elem = 18, num_loc = 1, k = 2;
int num_node = mesh.GetNodalNumber(num_elem, num_loc, k);

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetNormale

Syntax

void GetNormale(const Matrix2_2& dfjm1, R2& normale, Real_wp& dsj, int num_elem, int num_loc)

This method computes the normale and lineic element dsi from the inverse of jacobian matrix.

Parameters

dfjm1 (in)
inverse of jacobian matrix
normale (out)
outgoing normale to the selected edge
dsj (out)
lineic element dsi (used for integration)
num_elem (in)
element number
num_loc (in)
local edge number (in the element)

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

// then computation of dfj = DFi(pt_loc)
R2 pt_loc(1.0, 0.4);
Matrix2_2 dfjm1;
mesh.DFj(s, Pts, pt_loc, dfjm1, mesh, n);
GetInverse(dfjm1); // we store dfi^-1

// normale and ds
R2 normale; Real_wp ds;
int num_loc = 1; // local number of edge
mesh.GetNormale(dfjm1, normale, ds, n, num_loc);

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetDistanceToBoundary

Syntax

Real_wp GetDistanceToBoundary(const R2& pt_loc, int n)

This method returns the distance between a point and the boundary on the reference element. If the distance is negative, it means that the point is outside the reference element.

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

R2 pt_loc(0.1, 0.3); int n = 10;
Real_wp dist = mesh.GetDistanceToBoundary(pt_loc, n);

Location :

MeshBase.hxx
MeshBase.cxx

ProjectPointOnBoundary

Syntax

void ProjectPointOnBoundary(R2& pt_loc, int n)

This method projects a point on the boundary of reference element.

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

R2 pt_loc(0.01, 1.0);

// if the point is slightly outside the reference element, you may want to project it on the boundary
mesh.ProjectPointOnBoundary(pt_loc, n);

Location :

MeshBase.hxx
MeshBase.cxx

GetBoundingBox

Syntax

void GetBoundingBox(const VectR2& s, Real_wp coef, VectR2& box, TinyVector<R2, 2>& enveloppe)
void GetBoundingBox(int i, const VectR2& s, const SetPoints<Dimension2>& pts, TinyVector<R2, 2>& enveloppe)

This method computes the bounding box containing the element defined by its vertices s. box is a bounding polygon, enveloppe the bounding box, and coef is a safety coefficient (which can be used for curved elements for example). In the second syntax, the user provides the nodal points of the element, such that the computed bounding box is much more accurate.

Parameters

i (in)
element number
s (in)
vertices of the element
Pts (in)
nodal points
coef (in)
safety coefficient (for curved elements)
box (out)
polygonal enveloppe
enveloppe (out)
rectangular bounding box

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// bounding box of element 18 ?
VectR2 s;
int n = 18;
mesh.GetVerticesElement(n, s);

// estimating the bounding box with a safety coefficient of 10%
// put 0 if the element is straight
VectR2 box(s.GetM());
TinyVector<R2, 2> enveloppe;
mesh.GetBoundingBox(s, 0.1, box, enveloppe);
cout << "the element is included in box " << "[" << enveloppe(0)(0) << "," << enveloppe(1)(0) << "] x ["
     << enveloppe(0)(1) << "," << enveloppe(1)(1) << "]" << endl;

// other approach : compute nodal points
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

// in that case the bounding box is more accurate
mesh.GetBoundingBox(n, s, Pts, enveloppe);

Location :

MeshBase.hxx
MeshBase.cxx

ComputeValuesPhiNodalRef

Syntax

void ComputeValuesPhiNodalRef(int i, const R2& pt_loc, VectReal_wp& phi) const

This method computes the nodal shape functions on the local point pt_loc.

Parameters

i (in)
element number
pt_loc (in)
local coordinates of point
phi (out)
values of shape functions on the given point

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// you can compute the values of shape functions of element 18 on a local point
R2 pt_loc(0.6, 0.4); // local coordinates on the reference element
mesh.ComputeValuesPhiNodalRef(18, pt_loc, phi);

Location :

MeshBase.hxx
MeshBase.cxx

ComputeNodalGradient

Syntax

void ComputeNodalGradient(const SetMatrices<Dimension2>& mat, const Vector& Uloc, Vector& gradU, int i) const

This method computes the gradient of Uloc on nodal points of the element i.

Parameters

mat (in)
jacobian matrices on nodal points
Uloc (in)
values of u on nodal points
gradU (out)
gradient of u on nodal points
i (in)
element number

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

// nodal points of element 18
int n = 18;
VectR2 s;
mesh.GetVerticesElement(n, s);
SetPoints<Dimension2> Pts;
mesh.FjElemNodal(s, Pts, mesh, n);

SetMatrices<Dimension2> Mat;
mesh.DFjElemNodal(s, Pts, Mat, mesh, n);

// then if you define Uloc on the element
int nb_nodes = Pts.GetNbPointsNodal();
Vector<double> Uloc(nb_nodes);
double x, y;
for (int i = 0; i < nb_nodes; i++)
  {
    x = Pts.GetPointNodal(i)(0);
    y = Pts.GetPointNodal(i)(1);
    Uloc(i) = sin(pi_wp*x)*sin(pi_wp*y)*(x+y);
  }

// you can compute the gradient on nodal points of the element
VectR2 gradU(nb_nodes);
mesh.ComputeNodalGradient(Mat, Uloc, gradU, n);

Location :

MeshBase.hxx
MeshBase.cxx

GetEdgesOnZaxis

Syntax

void GetEdgesOnZaxis(IVect& list_vertices, IVect& list_edges, VectBool& Vertex_On_Axe)

This method finds all the edges and vertices on axis x = 0. It is useful for axisymmetric computation. list_vertices is an array that will contain all the vertex numbers of vertices located on this axis, list_edges contains edge numbers. Vertex_On_Axe(i) is equal to true if the vertex i is located on the axis.

Parameters

list_vertices (out)
numbers of vertices that are located on axis x = 0
list_edges (out)
numbers of edges that are located on axis x = 0
Vertex_On_Axe (out)
if Vertex_On_Axe(i) is true, the vertex is on the axis x = 0

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

IVect list_vertices, list_edges;
VectBool VertexOnAxe;
mesh.GetEdgesOnZaxis(list_vertices, list_edges, Vertex_On_Axe);

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

GetTetaMinMax

Syntax

void GetTetaMinMax(R2 center, Real_wp& teta_min, Real_wp& teta_max)

This method is used to determine the minimal and maximal value of theta (in polar coordinates) of all the vertices of the mesh.

Parameters

center (in)
center of polar coordinates
teta_min (out)
minimal value of angle θ
teta_max (out)
maximal value of angle θ

Example :

Mesh<Dimension2> mesh, mesh2;
mesh.Read("toto.mesh");

double teta_min, teta_max;
R2 center;
mesh.GetTetaMinMax(center, teta_min, teta_max);

Related topics :

AppendMesh

Location :

Mesh2D.hxx
Mesh2D.cxx

GetBoundaryRotation

Syntax

int GetBoundaryRotation(int f, int e1, int e2, int pos1, int& pos2, int& rot) const

When an edge f is shared by two elements e1 and e2, and we know that the local position of edge f is pos, this method finds the local position on the edge for the second element (pos2), and the difference of orientation of the edge between the two elements (rot). When the edge is located on a periodic boundary, it returns the edge number of the corresponding boundary.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

int num_edge = 10;
int num_elem1 = mesh.Boundary(num_edge).numElement(0);
int num_elem2 = mesh.Boundary(num_edge).numElement(1);
int pos1 = mesh.Element(num_elem1).GetPositionBoundary(num_edge);
int pos2, rot;
mesh_num.GetBoundaryRotation(num_edge, num_elem1, num_elem2, pos1, pos2, rot);

Location :

NumberMesh.hxx
NumberMesh.cxx

GetOrder

Syntax

int GetOrder() const
void GetOrder(TinyVector<IVect, 4>& order) const

This method returns the order of approximation. If a variable order is specified, it should return the maximal order of approximation. If you use the second syntax, you will get all the orders present in the mesh : order(0) will contain the orders for triangles, order(1) for quadrilaterals.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(3);
mesh_num.number_map.SetOrder(mesh, 3);


cout << "Order of approximation : " << mesh_num.GetOrder() << endl;

// for variable order, orders are retrieved in arrays
TinyVector<IVect, 4> order;
mesh_num.GetOrder(order);
// orders for triangles in order(0)
cout << "Order for triangles " << order(0) << endl;
// orders for quadrangles in order(1)
cout << "Order for quadrangles " << order(1) << endl;

Related topics :

SetOrder

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetNbDof

Syntax

int GetNbDof() const

This method returns the number of degrees of freedom contained in the mesh. If NumberMesh has not been called before, there will be no degree of freedom.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.NumberMesh();

cout << "Number of degrees of freedom : " << mesh_num.GetNbDof() << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetNbDof

Syntax

void SetNbDof(int N)

This method sets the number of degrees of freedom contained in the mesh. In a regular use, this number is set when NumberMesh is called and do not need to be modified.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.NumberMesh();

cout << "Number of degrees of freedom : " << mesh_num.GetNbDof() << endl;

// if you want to increase this value
int N = mesh_num.GetNbDof();
mesh_num.SetNbDof(N+1);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetNbLocalDof

Syntax

int GetNbLocalDof(int i) const

This method returns the number of degrees of freedom within an element of the mesh. If NumberMesh has not been called before, it should return 0.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.NumberMesh();

for (int i = 0; i < mesh.GetNbElt(); i++)
  cout << "Element i contains " << mesh_num.GetNbLocalDof(i) << "degrees of freedom " << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetNbDofPML

Syntax

int GetNbDofPML() const

This method returns the number of degrees of freedom contained in the PMLs. If NumberMesh has not been called before or if there are no pmls in the mesh, it will return 0.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

// you can add a PML
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0);

int nb_layers = 2; Real_wp delta = 0.5;
pml.SetAdditionPML(pml.PML_BOTH_SIDES, pml.PML_BOTH_SIDES, pml.PML_NO, nb_layers);
pml.SetThicknessPML(delta);
pml.AddPML(0, mesh);

// number the mesh
MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.compute_dof_pml = true;
mesh_num.NumberMesh();

// number of degrees of freedom for all the mesh :
cout << "Number of degrees of freedom : " << mesh_num.GetNbDof() << endl;

// number of degrees of freedom for the PML only :
cout << "Number of degrees of freedom in PML : " << mesh_num.GetNbDofPML() << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetDofPML

Syntax

int GetDofPML(int i) const

This method returns the PML dof number of the i-th dof. If the i-th dof is not in a PML, the method will return -1.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

// you can add a PML
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0);

int nb_layers = 2; Real_wp delta = 0.5;
pml.SetAdditionPML(pml.PML_BOTH_SIDES, pml.PML_BOTH_SIDES, pml.PML_NO, nb_layers);
pml.SetThicknessPML(delta);
pml.AddPML(0, mesh);

// number the mesh
MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.compute_dof_pml = true;
mesh_num.NumberMesh();

// number of degrees of freedom for all the mesh :
cout << "Number of degrees of freedom : " << mesh_num.GetNbDof() << endl;

// number of degrees of freedom for the PML only :
cout << "Number of degrees of freedom in PML : " << mesh_num.GetNbDofPML() << endl;

// loop over dofs to check dofs in PML
for (int i = 0; i < mesh_num.GetNbDof(); i++)
  if (mesh_num.GetDofPML(i) >= 0)
     cout << "Dof " << i << " is inside a PML region" << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ReallocateDofPML

Syntax

void ReallocateDofPML(int N)

This method changes the number of degrees of freedom in PML. Usually this number is already computed when the mesh is numbered and does not need to be changed.

Related topics :

GetNbDofPML

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetDofPML

Syntax

void SetDofPML(int i, int num_pml)

This method changes the dof PML number of the i-th dof. Usually this number is already computed when the mesh is numbered and does not need to be changed.

Related topics :

GetDofPML

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrder

Syntax

void SetOrder(int r)

This method sets the order of approximation associated with a mesh numbering. The same order will be set for all the elements of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

// setting an uniform order for all elements
MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);
mesh_num.NumberMesh();

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

IsOrderVariable

Syntax

bool IsOrderVariable() const

This method returns true if a variable order is set for the mesh, false otherwise. Variable order will be present if you have called ComputeVariableOrder or SetOrderElement.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");  

// setting an uniform order
MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// should return false
cout << "Mesh with variable order ? " << mesh_num.IsOrderVariable() << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetVariableOrder

Syntax

void SetVariableOrder(int type)

This method sets the type of algorithm used to attribute order to each edge, face, element of the mesh. Currently, there are two provided algorithms MEAN_EDGE_ORDER, and MAX_EDGE_ORDER. The first one uses the average length value of the edges to determine the size h of an element, whereas the second one takes the maximum edge length. The user can also provide the orders for each edge/element, in that case USER_ORDER must be chosen instead of MAX_EDGE_ORDER.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// type of variable order
// if you choose mesh_num.USER_ORDER, you have to provide the orders for each edge and face of the mesh
mesh_num.SetVariableOrder(mesh_num.MAX_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// setting the number of dofs for the orders computed
TinyVector<IVect, 4> order_used;
mesh_num.GetOrder(order_used);
for (int i = 0; i < order_used(0).GetM(); i++)
   mesh_num.number_map.SetOrder(mesh, order_used(0)(i));

// and numbering
mesh_num.NumberMesh();

// the orders are written on a mesh
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetMeshSizeVariableOrder

Syntax

void SetMeshSizeVariableOrder(const VectReal_wp& step)

This method sets the mesh sizes h(r) (given in the parameter step) used to attribute order to each edge. If the length of an edge is between h(r) and h(r+1), the order associated with this edge will be equal to r.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MAX_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// setting the number of dofs for the orders computed
TinyVector<IVect, 4> order_used;
mesh_num.GetOrder(order_used);
for (int i = 0; i < order_used(0).GetM(); i++)
   mesh_num.number_map.SetOrder(mesh, order_used(0)(i));

// and numbering
mesh_num.NumberMesh();

// the orders are written on a mesh
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetCoefficientVariableOrder

Syntax

void SetCoefficientSizeVariableOrder(const VectReal_wp& coef)

This method sets the coefficient attributed to each element in order to determine the order of each edge, face, element of the mesh. If the length of an edge multiplied by this coefficient is between h(r) and h(r+1), the order associated with this edge will be equal to r.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// then for each element you can attribute a coefficient
// usually we take into account velocities on each physical domain
Vector<double> coef(mesh.GetNbElt());
coef.Fill(1);
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.Element(i).GetReference() == 2)
    coef(i) = 2.0;

mesh_num.SetCoefficientVariableOrder(coef);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// setting the number of dofs for the orders computed
TinyVector<IVect, 4> order_used;
mesh_num.GetOrder(order_used);
for (int i = 0; i < order_used(0).GetM(); i++)
   mesh_num.number_map.SetOrder(mesh, order_used(0)(i));

// and numbering
mesh_num.NumberMesh();

// the orders are written on a mesh
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOrderElement

Syntax

int GetOrderElement(int i) const

This method returns the order of approximation for an element of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh.SetMeshSizeVariableOrder(step);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// printing orders of element ?
for (int i = 0; i < mesh.GetNbElt(); i++)
  cout << "Order for element " << i << " : " << mesh_num.GetOrderElement(i) << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOrderEdge

Syntax

int GetOrderEdge(int i) const

This method returns the order of approximation for an edge of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// printing orders of edges ?
for (int i = 0; i < mesh.GetNbEdges(); i++)
  cout << "Order for edge " << i << " : " << mesh_num.GetOrderEdge(i) << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOrderQuadrature

Syntax

int GetOrderQuadrature(int i) const
void GetOrderQuadrature(TinyVector<IVect, 4>& order) const

This method returns the order of quadrature formula for an edge of the mesh. It is usually set as the maximum order of the two adjoining elements so that boundary integrals are exactly evaluated.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(false);

// printing orders of quadrature ?
for (int i = 0; i < mesh.GetNbEdges(); i++)
  cout << "Order of quadrature for edge " << i << " : " << mesh.GetOrderQuadrature(i) << endl;

// if you want to know all the orders present in the mesh
TinyVector<IVect, 4> order_used;
mesh_num.GetOrderQuadrature(order_used);
cout << "Orders of quadrature present :  " << order_used << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOrderFace

Syntax

int GetOrderFace(int i) const

This method should not be used for 2-D meshes but only 3-D meshes.

GetOrderInside

Syntax

int GetOrderInside(int i) const

This method returns the order of approximation for the interior of element i of the mesh.

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrderElement

Syntax

void SetOrderElement(int i, int r)

This method sets the order of approximation for an element of the mesh. If the mesh was previously specified with a constant order, it will be switched to a variable order, so that you can manually attribute an order to each element of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// attributing orders to each element
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    int r = 3;
    if (i%2 == 1)
      r = 4;
   
    mesh_num.SetOrderElement(i, r);
  }

// and on edges
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  {
    int r = 1;
    for (int k = 0; k < mesh.Boundary(i).GetNbElements(); k++)
      {
        int ne = mesh.Boundary(i).numElement(k);
        r = max(r, mesh.GetOrderElement(ne));
      }
   
    mesh_num.SetOrderEdge(i, r);
  }

// setting the number of dofs for the orders computed
mesh_num.number_map.SetOrder(mesh, 3);
mesh_num.number_map.SetOrder(mesh, 4);

// and numbering
mesh_num.NumberMesh();

// the orders are written on a mesh
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrderEdge

Syntax

void SetOrderEdge(int i, int r)

This method sets the order of approximation for an edge of the mesh. If the mesh was previously specified with a constant order, it will be switched to a variable order, so that you can manually attribute an order to each element and edge of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// attributing orders to each element
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    int r = 3;
    if (i%2 == 1)
      r = 4;
   
    mesh_num.SetOrderElement(i, r);
  }

// and on edges
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  {
    int r = 1;
    for (int k = 0; k < mesh.Boundary(i).GetNbElements(); k++)
      {
        int ne = mesh.Boundary(i).numElement(k);
        r = max(r, mesh.GetOrderElement(ne));
      }
   
    mesh_num.SetOrderEdge(i, r);
  }

// setting the number of dofs for the orders computed
mesh_num.number_map.SetOrder(mesh, 3);
mesh_num.number_map.SetOrder(mesh, 4);

// and numbering
mesh_num.NumberMesh();

// the orders are written on a mesh
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrderQuadrature

Syntax

void SetOrderQuadrature(int i, int r)

This method sets the order of quadrature for an edge of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// attributing orders to each element
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    int r = 3;
    if (i%2 == 1)
      r = 4;
   
    mesh_num.SetOrderElement(i, r);
  }

// and quadrature on edges
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  {
    int r = 1;
    for (int k = 0; k < mesh.Boundary(i).GetNbElements(); k++)
      {
        int ne = mesh.Boundary(i).numElement(k);
        r = max(r, mesh.GetOrderElement(ne));
      }

    mesh_num.SetOrderQuadrature(i, r);
  }

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrderFace

Syntax

void SetOrderFace(int i, int r)

This method should not be used for 2-D meshes but only 3-D meshes.

SetOrderInside

Syntax

void SetOrderInside(int i, int r)

This method sets the order of approximation for the interior of element i of the mesh.

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetNbPointsQuadratureBoundary

Syntax

int GetNbPointsQuadratureBoundary(int i) const

This method returns the number of quadrature points used on an edge of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
mesh_num.ComputeVariableOrder(true);

// and numbering (assuming that mesh_num.number_map has been filled correctly)
mesh_num.NumberMesh();

// then you can know the number of quadrature points on each edge
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  cout << "Edge " << i << " contains " << mesh_num.GetNbPointsQuadratureBoundary(i) << " quadrature points" << endl;

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ClearDofs

Syntax

void ClearDofs()

This method clears degrees of freedom present in the mesh. It cancels task performed by method NumberMesh.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh.Read("toto.mesh");
mesh_num.SetOrder(5);

// mesh is numbered
mesh_num.NumberMesh();

// then you remove dof numbers
mesh_num.ClearDofs();

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ComputeVariableOrder

Syntax

void ComputeVariableOrder(int dg_formulation)

This method attributes an order of approximation to each edge and element of the mesh in the case where the order is variable. This should be called before NumberMesh so that numbering of the mesh is correctly performed. In 2-D, we attribute an order re to each edge e of the mesh with the following rule :

where he is the length of the edge e, steps λ are specified by method SetMeshSizeVariableOrder, and c_e is a coefficient controlled by method SetCoefficientVariableOrder. When orders have been attributed to the edges, we then attribute an order to each element of the mesh. If MEAN_EDGE_ORDER algorithm has been selected, the order for each element will be equal to

If MAX_EDGE_ORDER algorithm has been selected, the order for each element will be equal to

If dg_formulation is set to DISCONTINUOUS or HDG, the order of edges will then be modified to be the maximal value of orders of neighboring elements, so that quadrature formulas on edges of the mesh are always accurate enough. This is important for Discontinuous Galerkin formulation where boundary integrals defined on internal edges are present.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// defining the values for h (in step)
Vector<double> step(10);
step.Fill(0);

// setting intervals for each order
// if h < 0.05 => first order
step(2) = 0.05;
// if 0.05 <= h < 0.15 => second order
step(3) = 0.15;
// if 0.15 <= h < 0.3 => third order
step(4) = 0.3;
// etc
mesh_num.SetMeshSizeVariableOrder(step);

// type of variable order
mesh_num.SetVariableOrder(mesh_num.MEAN_EDGE_ORDER);

// then computation of order for each edge, face, element
// choose between CONTINUOUS, DISCONTINUOUS or HDG for the parameter
mesh_num.ComputeVariableOrder(ElementReference_Base::DISCONTINUOUS);

// and numbering (assuming that mesh_num.number_map is already filled)
mesh_num.NumberMesh();

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMesh.cxx

NumberMesh

Syntax

void NumberMesh()

This methods numbers the mesh, by attributing global dof numbers to each element. The number of dofs per vertex, edge, face and element is defined in object number_map. If discontinuous formulation is set in object number_map, we store only the first dof number on each element (since other dof numbers of the element are obtained by incrementing this initial dof number). Local numbers for reference elements are assumed to respect the following order :

The Montjoie finite element classes are based on that ordering.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);

// for specific finite element classes you can use ConstructNumberMap :
QuadrangleHcurlFirstFamily fe_hcurl;
fe_hcurl.ConstructFiniteElement(order);
fe_hcurl.ConstructNumberMap(mesh_num.number_map, false);

// the mesh is numbered
mesh_num.NumberMesh();

cout << "Number of degrees of freedom : " << mesh_num.GetNbDof() << endl;

Location :

NumberMesh.hxx
NumberMesh.cxx

GetNbDofElement

Syntax

int GetNbDofElement(int i, Edge<Dimension2> edge)
int GetNbDofElement(int i, Face<Dimension2> elt)

This method is used by NumberMesh to know the expected number of degrees of freedom for an element (or an edge). On a regular use, you should call GetNbLocalDof if you want to know the number of degrees of freedom for an element of the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("toto.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
int order = 3;
mesh_num.SetOrder(order);

// using classical nodal numbering scheme (one dof per vertex, r dofs per edge, etc)
mesh_num.number_map.SetOrder(mesh, order);

// for specific finite element classes you can use ConstructNumberMap :
QuadrangleHcurlFirstFamily fe_hcurl;
fe_hcurl.ConstructFiniteElement(order);
fe_hcurl.ConstructNumberMap(mesh_num.number_map, false);

// the mesh is numbered
mesh_num.NumberMesh();

// number of dofs on an edge of the mesh :
int num_edge = 42;
int nb_dof_edge = mesh_num.GetNbDofElement(num_edge, mesh.Boundary(num_edge));

// number of dofs on an element
int num_elt = 31;
int nb_dof_elt = mesh_num.GetNbDofElement(num_elt, mesh.Element(num_elt));

// for this last case, GetNbLocalDof is fastest :
nb_dof_elt = mesh_num.GetNbLocalDof(num_elt);

Location :

NumberMesh.hxx
NumberMesh.cxx

ConstructQuadrilateralNumbering

Syntax

void ConstructQuadrilateralNumbering(int r, Matrix<int>& NumNodes, Matrix<int>& CoordinateNodes)

For standard nodal quadrilateral finite element, this method gives you local numbers of the square as they should be ordered :

Loading

NumNodes and CoordinatesNodes allows you to switch from a 2-index notation (i1, i2) to a 1-index notation i, as shown in the example below :

Example :

int r = 3;
Matrix<int> NumNodes, CoordinateNodes;
MeshNumbering<Dimension2>::ConstructQuadrilateralNumbering(r, NumNodes, CoordinateNodes);

// you may want to loop over x-coordinate and y-coordinate
for (int i1 = 0; i1 <= r; i1++)
  for (int i2 = 0; i2 <= r; i2++)
    {
      // then retrieve the number associated with point (xi_i1, xi_i2)
      int node = NumNodes(i1, i2);
    }

// you may want to loop over nodes
for (int i = 0; i < (r+1)*(r+1); i++)
  {
    // and retrieve x-coordinate, y-coordinate of this node
    int i1 = CoordinateNodes(i, 0);
    int i2 = CoordinateNodes(i, 1);
  }

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMesh.cxx

ConstructTriangularNumbering

Syntax

void ConstructTriangularNumbering(int r, Matrix<int>& NumNodes, Matrix<int>& CoordinateNodes)

For a standard nodal triangular finite element, this method gives you local numbers of the triangle as they should be ordered :

Loading

NumNodes and CoordinatesNodes allows you to switch from a 2-index notation (i1, i2) to a 1-index notation i, as shown in the example below :

Example :

int r = 4;
Matrix<int> NumNodes, CoordinateNodes;
MeshNumbering<Dimension2>::ConstructTriangularNumbering(r, NumNodes, CoordinateNodes);

// you may want to loop over x-coordinate and y-coordinate
for (int i1 = 0; i1 <= r; i1++)
  for (int i2 = 0; i2 <= r-i1; i2++)
    {
      // then retrieve the number associated with point (xi_i1, xi_i2)
      int node = NumNodes(i1, i2);
    }

// you may want to loop over nodes
for (int i = 0; i < (r+1)*(r+2)/2; i++)
  {
    // and retrieve x-coordinate, y-coordinate of this node
    int i1 = CoordinateNodes(i, 0);
    int i2 = CoordinateNodes(i, 1);
  }

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMesh.cxx

GetRotationFace

Syntax

int GetRotationFace(int way1, int way2, int nv)

This method will return 0 if way1 and way2 are equal (edges with same orientation), and 1 if way1 and way2 are different (edges with opposite orientation).

Parameters

way1
orientation of first edge
way2
orientation of second edge
nv
number of vertices (unused)

Example :

MeshNumbering<Dimension2> mesh_num(mesh);
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  if (mesh.Boundary(i).GetNbElements() == 2)
    {
      int n1 = mesh.Boundary(i).numElement(0);
      int n2 = mesh.Boundary(i).numElement(1);
      int pos1 = mesh.Element(n1).GetPositionBoundary(i);
      int pos2 = mesh.Element(n2).GetPositionBoundary(i);
      int way1 = mesh.Element(n1).GetOrientationBoundary(pos1);
      int way2 = mesh.Element(n2).GetOrientationBoundary(pos2);
      int rot = mesh_num.GetRotationFace(way1, way2, 2);
    }

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMesh.cxx

GetOppositeOrientationFace

Syntax

int GetOppositeRotationFace(int rot, Edge<Dimension2>& e)

This method will return rot.

Example :

MeshNumbering<Dimension2> mesh_num(mesh);
for (int i = 0; i < mesh.GetNbBoundary(); i++)
  if (mesh.Boundary(i).GetNbElements() == 2)
    {
      int n1 = mesh.Boundary(i).numElement(0);
      int n2 = mesh.Boundary(i).numElement(1);
      int pos1 = mesh.Element(n1).GetPositionBoundary(i);
      int pos2 = mesh.Element(n2).GetPositionBoundary(i);
      int way1 = mesh.Element(n1).GetOrientationBoundary(pos1);
      int way2 = mesh.Element(n2).GetOrientationBoundary(pos2);
      int rot = mesh.GetRotationFace(way1, way2, 2);
      // opposite orientation (but in 2-D, it is the same)
      int rot2 = mesh_num.GetOppositeOrientationFace(rot, mesh.Boundary(i));
    }

Related topics :

NumberMesh

Location :

NumberMesh.hxx
NumberMesh.cxx

GetNbReferences

Syntax

int GetNbReferences() const

This method returns the number of references (Physical Line for Gmsh) potentially present in the mesh. We take into account only references related to edges. If you want to create a new reference number, use the method GetNewReference so that arrays are resized correctly.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// loop over references of the mesh
for (int ref = 0; ref <= mesh.GetNbReferences(); ref++)
  if (mesh.GetCurveType(ref) > 0)
    cout << "Type of curve : " << mesh.GetCurveType(ref) << endl;

Related topics :

GetNewReference
ResizeNbReferences

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetCurveType

Syntax

int GetCurveType(int ref) const

This method returns the curve type associated with a reference (e.g. CURVE_CIRCLE, CURVE_SPLINE).

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// loop over references of the mesh
for (int ref = 0; ref <= mesh.GetNbReferences(); ref++)
  if (mesh.GetCurveType(ref) > 0)
    cout << "Type of curve : " << mesh.GetCurveType(ref) << endl;

Related topics :

SetCurveType

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetBoundaryCondition

Syntax

int GetBoundaryCondition(int ref) const
const IVect& GetBoundaryCondition() const

This method returns the boundary condition associated with a reference (e.g. LINE_DIRICHLET, LINE_NEUMANN).

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// loop over references of the mesh
for (int ref = 0; ref <= mesh.GetNbReferences(); ref++)
  if (mesh.GetBoundaryCondition(ref) > 0)
    cout << "Boundary condition of curve : " << mesh.GetBoundaryCondition(ref) << endl;

// you can also retrieve the array containing all the boundary conditions
const IVect& Cond_curve = mesh.GetBoundaryCondition();

Related topics :

SetBoundaryCondition

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetBodyNumber

Syntax

int GetBodyNumber(int ref) const
const IVect& GetBodyNumber() const

This method returns the body number associated with a reference. Bodies are currently used to specify transparent condition in Montjoie.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// you can get a single body number
int n = mesh.GetBodyNumber(4);

// you can also retrieve the array containing all the body numbers
const IVect& Cond_curve = mesh.GetBodyNumber();

Related topics :

SetBodyNumber

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetCurveParameter

Syntax

void GetCurveParameter(int ref, VectReal_wp& param) const
const VectReal_wp& GetCurveParameter(int ref) const

This method returns the curve parameters associated with a curve. For example, parameters of circles consist of the center and radius.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
mesh.Read("test.mesh");

// the parameters of circle should have been found 
Vector<double> param = mesh.GetCurveParameter(1);

// center of circle and radius
R2 center; Real_wp radius;
center.Init(param(0), param(1));
radius = param(2); 

Related topics :

SetCurveParameter

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetCurveFunctions1D

Syntax

Globatto<Real_wp>& GetCurveFunctions1D() const

This method returns the basis functions associated with Gauss-Lobatto points (for curved edges).

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetOrderGeometry(4);
mesh.Read("test.msh");

// if you need basis functions used for curved boundaries
Globatto<Real_wp>& lob = mesh.GetCurveFunctions1D();

// you can evaluate basis functions
VectReal_wp phi; Real_wp x(0.6);
lob.ComputeValuesPhiRef(x, phi);
DISP(phi);

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetNbPmlAreas

Syntax

int GetNbPmlAreas() const

This method returns the number of PML regions in the mesh.

Example :

Mesh<Dimension2> mesh;

// you change the number of PML regions
mesh.ReallocatePmlAreas();

// you can access to a PML region
PmlRegionParameter<Dimension2>& pml0 = mesh.GetPmlArea(0);
PmlRegionParameter<Dimension2>& pml1 = mesh.GetPmlArea(1);

// sets PML parameters for each region
pml0.SetAdditionPML(pml0.PML_POSITIVE_SIDE, pml0.PML_NO, pml0.PML_NO, 3);
pml1.SetAdditionPML(pml0.PML_NO, pml0.PML_NEGATIVE_SIDE, pml0.PML_NO, 4);

// if you call ConstructMesh, the PMLs will be automatically added
// with the provided parameters
Vector<string> param(1); param(0) = "test.mesh";
mesh.ConstructMesh(mesh.QUADRILATERAL_MESH, param);

// you can retrieve the number of PML regions
int nb_pml_area = mesh.GetNbPmlAreas()

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

ReallocatePmlAreas

Syntax

void ReallocatePmlAreas(int n)

This method changes the number of PML regions in the mesh (previous parameters are cleared).

Example :

Mesh<Dimension2> mesh;

// you change the number of PML regions
mesh.ReallocatePmlAreas();

// you can access to a PML region
PmlRegionParameter<Dimension2>& pml0 = mesh.GetPmlArea(0);
PmlRegionParameter<Dimension2>& pml1 = mesh.GetPmlArea(1);

// sets PML parameters for each region
pml0.SetAdditionPML(pml0.PML_POSITIVE_SIDE, pml0.PML_NO, pml0.PML_NO, 3);
pml1.SetAdditionPML(pml0.PML_NO, pml0.PML_NEGATIVE_SIDE, pml0.PML_NO, 4);

// if you call ConstructMesh, the PMLs will be automatically added
// with the provided parameters
Vector<string> param(1); param(0) = "test.mesh";
mesh.ConstructMesh(mesh.QUADRILATERAL_MESH, param);

// you can retrieve the number of PML regions
int nb_pml_area = mesh.GetNbPmlAreas()

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetPmlArea

Syntax

PmlRegionParameter<Dimension2>& GetPmlArea(int i) const

This method gives an access to a PML region.

Example :

Mesh<Dimension2> mesh;

// you change the number of PML regions
mesh.ReallocatePmlAreas();

// you can access to a PML region
PmlRegionParameter<Dimension2>& pml0 = mesh.GetPmlArea(0);
PmlRegionParameter<Dimension2>& pml1 = mesh.GetPmlArea(1);

// sets PML parameters for each region
pml0.SetAdditionPML(pml0.PML_POSITIVE_SIDE, pml0.PML_NO, pml0.PML_NO, 3);
pml1.SetAdditionPML(pml0.PML_NO, pml0.PML_NEGATIVE_SIDE, pml0.PML_NO, 4);

// if you call ConstructMesh, the PMLs will be automatically added
// with the provided parameters
Vector<string> param(1); param(0) = "test.mesh";
mesh.ConstructMesh(mesh.QUADRILATERAL_MESH, param);

// you can retrieve the number of PML regions
int nb_pml_area = mesh.GetNbPmlAreas()

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetBoundaryCondition

Syntax

void SetBoundaryCondition(int ref, int cond)

This method sets the boundary condition associated with a reference (e.g. LINE_DIRICHLET, LINE_NEUMANN). The list of boundary conditions is present in class BoundaryConditionEnum, the following boundary conditions are available :

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// setting Dirichlet boundary condition with reference 2
mesh.SetBoundaryCondition(2, mesh.LINE_DIRICHLET);

Related topics :

GetBoundaryCondition

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetCurveType

Syntax

void SetCurveType(int ref, int type)

This method sets the curve type associated with a reference (e.g. CURVE_CIRCLE, CURVE_SPLINE). The following types of curve are available

For CURVE_FILE, it is automatically managed (you don't need to specify it) when you read a file. For all the curves, it is required to call SetGeometryOrder with the wished order of approximation (for computation of intermediary points).

Example :

Mesh<Dimension2> mesh;

// setting order of geometry
mesh.SetGeometryOrder(3);

// setting curves before reading mesh
mesh.SetCurveType(2, mesh.CURVE_SPLINE);

// we read the mesh
mesh.Read("test.mesh");

Related topics :

GetCurveType

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetBodyNumber

Syntax

void SetBodyNumber(int ref, int num)

This method sets the body number associated with a reference. It may be useful to group several references into a single body. Bodies are currently used to specify transparent condition in Montjoie.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// you can set a body with references 2 4 5
int num = 1;
mesh.SetBodyNumber(2, num);
mesh.SetBodyNumber(4, num);
mesh.SetBodyNumber(5, num);

Related topics :

GetBodyNumber

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetCurveParameter

Syntax

void SetCurveParameter(int ref, const VectReal_wp& param)

This method modifies the parameters associated with a curve. For example, parameters of circles consist of the center and radius. Here the list of 2-D curves with parameters :

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
// parameters : center of circle and radius
Vector<double> param(3);
param(0) = 0;
param(1) = 0;
param(2) = 2.0;
mesh.SetCurveParameter(1, param);
mesh.Read("test.mesh");

Related topics :

GetCurveParameter

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

CopyCurves

Syntax

void CopyCurves(const Mesh<Dimension2>& mesh)

This method copies the curves of a mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
// parameters : center of circle and radius
Vector<double> param(3);
param(0) = 0;
param(1) = 0;
param(2) = 2.0;
mesh.SetCurveParameter(1, param);
mesh.Read("test.mesh");

Mesh<Dimension2> mesh2;

// we copy curves of first mesh
mesh2.CopyCurves(mesh);

Related topics :

GetCurveParameter

Location :

MeshBoundaries.hxx
MeshBoundaries.cxx

CopyInputData

Syntax

void CopyInputData(Mesh<Dimension2>& mesh)

This method copies the input datas of another mesh. The input datas are the parameters provided in the data file.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
// parameters : center of circle and radius
Vector<double> param(3);
param(0) = 0;
param(1) = 0;
param(2) = 2.0;
mesh.SetCurveParameter(1, param);
mesh.Read("test.mesh");

Mesh<Dimension2> mesh2;

// we copy all input data of the mesh
mesh2.CopyInputData(mesh);

Location :

MeshBoundaries.hxx
MeshBoundaries.cxx

GetNewReference

Syntax

int GetNewReference()

This method returns an unused reference number.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// new reference number
int new_ref = mesh.GetNewReference();

Related topics :

GetNbReferences

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

ResizeNbReferences

Syntax

void ResizeNbReferences(int n)

This method modifies the number of available references (arrays are resized). If you want to add a new reference, you can call GetNewReference instead.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

mesh.ResizeNbReferences(20);

Related topics :

GetNbReferences

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetNewReference

Syntax

void SetNewReference(int ref, int type, int cond)

This method modifies the curve type and boundary condition associated with a reference.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// new reference number
int new_ref = mesh.GetNewReference();

// setting curve type and boundary condition
mesh.SetNewReference(new_ref, mesh.CURVE_CIRCLE, mesh.LINE_NEUMANN);

Related topics :

GetNbReferences

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetThicknessPML

Syntax

Real_wp GetThicknessPML() const

This method returns the PML thickness.

Example :

PmlRegionParameter<Dimension2> pml;

Real_wp delta = pml.GetThicknessPML();

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetOriginRadialPML

Syntax

R2 GetOriginRadialPML() const

This method returns the center of PML.

Example :

PmlRegionParameter<Dimension2> pml;

R2 center = pml.GetOriginRadialPML();

Related topics :

SetOriginRadialPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetRadiusPML

Syntax

Real_wp GetRadiusPML() const

This method returns the PML radius.

Example :

PmlRegionParameter<Dimension2> pml;

Real_wp r = pml.GetRadiusPML();

Related topics :

SetRadiusPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetThicknessPML

Syntax

void SetThicknessPML(Real_wp delta)

This method sets PML thickness.

Example :

PmlRegionParameter<Dimension2> pml;

pml.SetThicknessPML(2.5);

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetTypeAdditionPML

Syntax

int GetTypeAdditionPML(int num_coor) const

This method returns PML type with respect to x-coordinate (or y-coordinate).

Example :

PmlRegionParameter<Dimension2> pml;

// pml added for x < xmin and/or x > xmax ?
int type_add_x = pml.GetTypeAdditionPML(0);

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetNbLayersPML

Syntax

int GetNbLayersPML() const

This method returns the numbers of layers in PML.

Example :

PmlRegionParameter<Dimension2> pml;

// number of layers ?
int n = pml.GetNbLayersPML();

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetRadiusPML

Syntax

void SetRadiusPML(Real_wp r)

This method sets the radius of PML.

Example :

PmlRegionParameter<Dimension2> pml;

// radius for PML layers
pml.SetRadiusPML(2.0);

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetOriginRadialPML

Syntax

void SetOriginRadialPML(R2 center)

This method sets the center of PML.

Example :

PmlRegionParameter<Dimension2> pml;

// center for PML layers
pml.SetOriginRadialPML(R2(1, 0));

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetAdditionPML

Syntax

void SetAdditionPML(int add_x, int add_y, int add_z, int nb_layers)

This method is used to set PML to add. PML are effectively added when you call AddPML. On x-coordinate, you can add PML on negative side (x < xmin), positive side (x > xmax) or both sides. If you set nb_layers to 0, the number of layers will be automatically computed when you call AddPML.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0);

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(pml.PML_NEGATIVE_SIDE, pml.PML_BOTH_SIDES, pml.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// mesh is modified in order to add those PML
pml.AddPML(0, mesh);

Related topics :

GetTypeAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundaries.cxx

SetParameters

Syntax

void SetParameters(VectString param)

This method sets the PML region with a line of the data file. This function is called when the data file is read (keyword AddPML). Outside of a data file, it seems easier to call SetAdditionPML.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0);

// adding PML for x < xmin, y < ymin and y > ymax
// here we use the syntax of AddPML keyword
Vector<string> param(4);
param(0) = "YES"; param(1) = "X-Y"; param(2) = to_str(2.0); param(3) = to_str(4);
pml.SetParameters(param);

// mesh is modified in order to add those PML
pml.AddPML(0, mesh);

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundaries.cxx

AddPML

Syntax

void AddPML(int num_region, Mesh<Dimension2>& mesh)

This method is used to add PML elements to a given mesh. num_region is the number of the PML region.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0); // number of the region is 0

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(pml.PML_NEGATIVE_SIDE, pml.PML_BOTH_SIDES, pml.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// mesh is modified in order to add those PML
pml.AddPML(0, mesh);

Related topics :

SetAdditionPML

Location :

MeshBoundaries.hxx
MeshBoundaries.cxx

GetNbPeriodicReferences

Syntax

int GetNbPeriodicReferences() const

This method returns the number of periodic conditions, each time you wrote PERIODICITY, PERIODIC_X, PERIODIC_Y, PERIODIC_Z or CYCLIC in the data file.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

// GetNbPeriodicReferences should return 1 here
int nb_per = mesh.GetNbPeriodicReferences();

Related topics :

AddPeriodicCondition

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetPeriodicReference

Syntax

int GetPeriodicReference(int num, int j) const

This method returns the first or second reference for the periodic condition num. Indeed, for each boundary condition, the user provides two references (the two curves that face them-selves), which can be accessed through this method.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, int n0 = mesh.GetPeriodicReference(0, 0);PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

// GetPeriodicReference to retrieve the numbers given in AddPeriodicCondition
int n0 = mesh.GetPeriodicReference(0, 0);
int n1 = mesh.GetPeriodicReference(0, 1);
// here n0 should be equal to 3 and n1 to 6

Related topics :

AddPeriodicCondition

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetSameNumberPeriodicDofs

Syntax

void SetSameNumberPeriodicDofs(bool same_number)

By calling this method, you can require (or not) that degrees of freedom located on periodic boundaries are attributed with the same numbers as the original degrees of freedom. If same numbers are attributed, there are no periodic dofs. The periodicity is ensured by construction, however quasi-periodic conditions are currently not implemented with this option. For quasi-periodic conditions (e.g. PERIODIC_X, PERIODIC_THETA), you cannot attribute same numbers.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
mesh_num.SetSameNumberPeriodicDofs(true);
mesh_num.NumberMesh();

Related topics :

UseSameNumberForPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

SetFormulationForPeriodicCondition, GetFormulationForPeriodicCondition

Syntax

void SetFormulationForPeriodicCondition(int type)
int GetFormulationForPeriodicCondition() const

The method SetFormulationForPeriodicCondition sets the formulation that will be used in Montjoie to handle periodic conditions (or quasi-periodic). You can choose between

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, PERIODIC_X);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
mesh_num.SetFormulationForPeriodicCondition(mesh_num.WEAK_PERIODIC);
mesh_num.NumberMesh();

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

UseSameNumberForPeriodicDofs

Syntax

bool UseSameNumberPeriodicDofs() const

This method returns true if degrees of freedom located on periodic boundaries are attributed with the same numbers as the original degrees of freedom.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(num);

// numbering the mesh (assuming that number_map is correctly filled)
cout << "Same number ? : " mesh_num.UseSameNumberForPeriodicDofs() << endl;
mesh_num.NumberMesh();

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

AddPeriodicCondition

Syntax

void AddPeriodicCondition(const TinyVector<int, 2>& ref, int type_periodicity)

This method adds periodic condition to the mesh. You give the two opposite references in the array ref, and the type of periodicity. The following types of periodicity are accepted :

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions for x-coordinate
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_X);

// and y-coordinate
num(0) = 2; num(1) = 5;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_Y);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);
// numbering the mesh (assuming that mesh_num.number_map is correctly filled)
mesh_num.NumberMesh();

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetNbPeriodicDof

Syntax

int GetNbPeriodicDof() const

This method returns the number of periodic dofs (obtained with a translation/rotation of original dofs).

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

cout << "Number of periodic dofs " << mesh_num.GetNbPeriodicDof() << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetPeriodicDof, SetPeriodicDof

Syntax

int GetPeriodicDof(int i) const
void SetPeriodicDof(int i, int num)

The method GetPeriodicDof returns the dof number of the i-th periodic dof of the mesh. The method SetPeriodicDof allows the user to change this number.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// display periodic and original dofs
for (int i = 0; i < mesh_num.GetNbPeriodicDof(); i++)
  cout << "Dof " << mesh_num.GetPeriodicDof(i) <<
       "obtained with a translation/rotation from dof " << mesh_num.GetOriginalPeriodicDof(i) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetPeriodicBoundary, SetPeriodicBoundary

Syntax

int GetPeriodicBoundary(int num_edge) const
void SetPeriodicBoundary(int n1, int n2)

The method GetPeriodicBoundary returns the edge number that is periodic with the given edge. If the edge does not belong to a periodic boundary, it should return -1. The method SetPeriodicBoundary allows the user to change this number.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// opposite edge to a given edge
int num_edge = 12;
int num_edge_opp = mesh_num.GetPeriodicBoundary(num_edge);
// if the edge num_edge has a reference 3, the edge num_edge_opp will have a reference 6 (ref 3 and 6 are said periodic)
// if num_edge_opp = -1, it means that the edge 12 did not belong to a periodic reference

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

InitPeriodicBoundary

Syntax

void InitPeriodicBoundary()

This method allocates arrays needed to store periodic boundaries. After calling this methode, you can call the method SetPeriodicBoundary to specify a periodic boundary. In a regular use, this method does not need to be called, since TreatPeriodicCondition fills the arrays for periodic boundaries.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// you did not number the mesh but you want to provide periodic boundaries
mesh_num.InitPeriodicBoundary();
mesh_num.SetPeriodicBoundary(3, 7);
// etc

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOriginalPeriodicDof, SetOriginalPeriodicDof

Syntax

int GetOriginalPeriodicDof(int i) const
void SetOriginalPeriodicDof(int i, int num)

The method GetOriginalPeriodicDof returns the original dof number of the i-th periodic dof of the mesh. The method SetOriginalPeriodicDof allows the user to change this number.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// display periodic and original dofs
for (int i = 0; i < mesh_num.GetNbPeriodicDof(); i++)
  cout << "Dof " << mesh_num.GetPeriodicDof(i) <<
       "obtained with a translation/rotation from dof " << mesh_num.GetOriginalPeriodicDof(i) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetTranslationPeriodicDof, SetTranslationPeriodicDof

Syntax

R2 GetTranslationPeriodicDof(int i) const
void SetTranslationPeriodicDof(int i, R2 u)

The method GetTranslationPeriodicDof returns the translation vector between two periodic degrees of freedom. The method SetTranslationPeriodicDof allows the user to change this translation vector.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// display translation
for (int i = 0; i < mesh_num.GetNbPeriodicDof(); i++)
  cout << "Dof " << mesh_num.GetPeriodicDof(i) <<
       "obtained with a translation of " << mesh_num.GetTranslationPeriodicDof(i) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetTranslationPeriodicBoundary

Syntax

R2 GetTranslationPeriodicBoundary(int i) const
VectR2& GetTranslationPeriodicBoundary() const

This method returns the translation vector between two periodic boundaries. With the second syntax, you can retrieve the translation vectors for all the referenced edges.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that mesh_num.number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// display translation vector between edge 0 and its opposite edge
cout << "Translation = " << mesh.GetTranslationPeriodicBoundary(0) << endl;

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetPeriodicityTypeReference

Syntax

int GetPeriodicityTypeReference(int i) const

This method returns the type of periodicity (in x, y, z or theta).

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// loop over periodic conditions
for (int i = 0; i < mesh.GetNbPeriodicReference(); i++)
  cout << "Type of periodicity for condition " << i << " : " << mesh.GetPeriodicityTypeReference(i) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetPeriodicityTypeForDof, SetPeriodicityTypeForDof

Syntax

int GetPeriodicityTypeForDof(int i) const
void SetPeriodicityTypeForDof(int i, int type)

The method GetPeriodicityTypeForDof returns the type of periodicity associated with the i-th periodic dof. It can be PERIODIC_CTE, PERIODIC_THETA, PERIODIC_X, PERIODIC_Y, PERIODIC_Z and also a combination of those : PERIODIC_XY, PERIODIC_XZ, PERIODIC_YZ, PERIODIC_XYZ, PERIODIC_ZTHETA. The method SetPeriodicityTypeForDof allows the user to change the type of periodicity.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
  
// adding periodic conditions
TinyVector<int, 2> num(2);
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_X);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// loop over periodic dofs
for (int i = 0; i < mesh_num.GetNbPeriodicDof(); i++)
  cout << "Type of periodicity for periodic dof " << i
       << " : " << mesh.GetPeriodicityTypeForDof(i) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetPeriodicityTypeForBoundary, SetPeriodicityTypeForBoundary

Syntax

int GetPeriodicityTypeForBoundary(int num_edge) const
void SetPeriodicityTypeForBoundary(int num_edge, int type)

The method GetPeriodicityTypeForBoundary returns the type of periodicity associated with the i-th edge. It can be PERIODIC_CTE, PERIODIC_THETA, PERIODIC_X, PERIODIC_Y, PERIODIC_Z. The method SetPeriodicityTypeForBoundary allows the user to change the type of periodicity.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
  
// adding periodic conditions
TinyVector<int, 2> num(2);
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_X);

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// numbering the mesh (assuming that number_map is correctly filled)
// if same dof numbers are used, no periodic dof
mesh_num.SetSameNumberPeriodicDofs(false);
mesh_num.NumberMesh();

// edge 12 is associated with a periodic condition ?
int type_per = mesh_num.GetPeriodicityTypeForBoundary(12);

Related topics :

SetSameNumberPeriodicDofs

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetPeriodicAlpha, SetPeriodicAlpha

Syntax

Real_wp GetPeriodicAlpha() const
void SetPeriodicAlpha(Real_wp alpha)

The method GetPeriodicAlpha returns the angle of the section (i.e. 2 Pi/N). The method SetPeriodicAlpha sets this angle. Usually, you don't need to call this method because the angle is automatically computed when the mesh is read (or constructed).

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions for cyclic mesh
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_THETA);

// reading the mesh
mesh.Read("exemple.mesh");

// when the mesh is read, the angle alpha should have been computed
double alpha = mesh.GetPeriodicAlpha();

Related topics :

AddPeriodicCondition

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

ClearPeriodicCondition

Syntax

void ClearPeriodicCondition()

This method clears informations about periodic conditions associated with the mesh.

Example :

Mesh<Dimension2> mesh;

// adding periodic conditions for cyclic mesh
TinyVector<int, 2> num;
num(0) = 3; num(1) = 6;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_THETA);

// you may want to remove them
mesh.ClearPeriodicConditions();

// and adding correct conditions
num(0) = 3; num(1) = 5;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// reading the mesh
mesh.Read("exemple.mesh");

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetSafetyCoef

Syntax

Real_wp GetSafetyCoef(int i) const

This method returns the safety coefficient used for an element in order to invert transformation Fi. Curved elements are attributed a non-null coefficient so that all the points of the display grid are correctly found.

Example :

Mesh<Dimension2> mesh;

cout << "Safety coefficient for element 0 = " << mesh.GetSafetyCoef(0) << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

AreCurveEqual

Syntax

bool AreCurveEqual(int ref1, int ref2) const

This method returns true if the curves associated with references ref1 and ref2 are the same.

Example :

Mesh<Dimension2> mesh;

int ref1 = 2;
int ref2 = 4;
mesh.SetCurveType(ref1, mesh.CURVE_CIRCLE);
mesh.SetCurveType(ref2, mesh.CURVE_CIRCLE);

mesh.Read("test.mesh");

if (mesh.AreCurveEqual(ref1, ref2)
  cout << "Curve equal" << endl;

Related topics :

SetSameNumberPeriodicDofs

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetGeometryOrder

Syntax

int GetGeometryOrder() const

This method returns the order of approximation r used for the definition of geometry. Each curved edge will then be described by two extremities and r-1 intermediary points (placed at Gauss-Lobatto points).

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// should give r = 3
int r = mesh.GetGeometryOrder();

Related topics :

SetGeometryOrder

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetInterpolate1D

Syntax

Real_wp GetInterpolate1D(int i, Real_wp x) const

This method returns an evaluation of Lagrange interpolation functions (based on Gauss-Lobatto points).

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// returns phi_0(x)
double x = 0.23;
double y = mesh.GetInterpolate1D(0, x);

Related topics :

SetGeometryOrder

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

GetNodalPointInsideEdge

Syntax

Real_wp GetNodalPointInsideEdge(int i) const

This method returns intermediary points used for curved elements. It should be Gauss-Lobatto points.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// intermediary points (on unit interval [0,1])
for (int i = 0; i < mesh.GetGeometryOrder()-1; i++)
  cout << i << "-th intermediary point " << mesh.GetNodalPointInsideEdge(i) << endl;

Related topics :

SetGeometryOrder

Location :

MeshBoundaries.hxx
MeshBoundariesInline.cxx

IsBoundaryCurved

Syntax

bool IsBoundaryCurved(int i) const

This method returns true if the i-th edge is curved.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.msh");

// edge 11 is curved ?
if (mesh.IsBoundaryCurved(11))
  cout << "Edge 11 is curved" << endl;

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

GetPointInsideEdge

Syntax

R2 GetPointInsideEdge(int i, int k) const

This method returns the k-th intermediary point of edge i.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// loop over referenced boundary
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  for (int k = 0; k < mesh.GetGeometryOrder()-1; k++)
    {
      // intermediary point of edge i
      R2 point = mesh.GetPointInsideEdge(i, k);
    }

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

SetPointInsideEdge

Syntax

void SetPointInsideEdge(int i, int k, R2 point)

This method sets the k-th intermediary point of edge i. This low-level routine should not be used in a regular use since, since these points are already computed when the mesh is constructed (or read from a file).

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// loop over referenced edges
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  {
    // extremities of the edge
    R2 ptA = mesh.Vertex(mesh.BoundaryRef(i).numVertex(0));
    R2 ptB = mesh.Vertex(mesh.BoundaryRef(i).numVertex(1));
    for (int k = 0; k < mesh.GetGeometryOrder()-1; k++)
      {
        // intermediary point of edge i, for example barycenter (straight edge)
        Real_wp lambda = mesh.GetNodalPointInsideEdge(k);
        R2 point = (1.0-lambda)*ptA + lambda*ptB;
        mesh.SetPointInsideEdge(i, k, point);
      }
   }

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

FindFollowingVertex

Syntax

int FindFollowingVertex(int n1, int n2, int ref, const IVect& ref_cond, int& e)
int FindFollowingVertex(int n1, int n2, int ref, const IVect& ref_cond, VectBool& edge_used, int& e)

This method finds the next vertex of the curve beginning with vertex n1 and n2. The curve can be a group of curves for which ref_cond(i) is equal to ref. The argument e is modified to contain the next edge number. You can also forbid some edges through array edge_used (if edge_used(i) is equal to true, edge i will not be considered when searching the next point).

Parameters

n1 (in)
first vertex of the curve
n2 (in)
second vertex of the curve
ref (in)
target number
ref_cond (in)
the curve is made of references i such that ref_cond(i) = ref
edge_used (in)
if edge_used(i) is true, the edge i is not considered
e (out)
edge number that links n2 to the next vertex

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// creating group of curves with reference 2 and 3
IVect ref_cond(mesh.GetNbReferences()+1);
ref_cond.Fill(0);
ref_cond(2) = 1;
ref_cond(3) = 1;

// finding an initial edge
int ne = -1;
for (int i = 0; i < mesh.GetNbBoundaryRef(); i++)
  if (ref_cond(mesh.BoundaryRef(i).GetReference()) == 1)
    ne = i;

// finding all points of the curve
int n1 = mesh.BoundaryRef(ne).numVertex(0);
int n2 = mesh.BoundaryRef(ne).numVertex(1);
int num_init = n1;
while (n2 != num_init)
  {
     // next vertex
     int inext = FindFollowingVertex(n1, n2, 1, ref_cond, ne);
     n1 = n2;
     n2 = inext;
   }

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

GetSurfaceFiniteElement

Syntax

EdgeLobatto& GetSurfaceFiniteElement()

This method returns the 1-D finite element used for curved boundaries.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.msh");

// we retrieve the 1-D finite element
EdgeLobatto& fe = mesh.GetSurfaceFiniteElement();

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

FindEdgeRefWithExtremities, FindEdgeWithExtremities

Syntax

int FindEdgeWithExtremities(int n1, int n2) const
int FindEdgeRefWithExtremities(int n1, int n2) const

This method returns the edge number of edge whose extremities are n1 and n2. You can also retrieve the reference edge number (but it is the same as global edge number, except if the edge is not referenced. In that case, the method will return -1). The cost of this method is in constant time, since we are looking all the edges leaving vertex n1 and n2.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// edge between vertex 13 and 18 ?
int ne = mesh.FindEdgeWithExtremities(13, 18);
if (ne < 0)
  cout << "Edge not existing " << endl;

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

FindVertexNumber

Syntax

int FindVertexNumber(const R2& pt)

This method returns the vertex number by knowing only its coordinates. The cost of this method is in linear time, since we are browsing all the vertices of the mesh.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// vertex number of point -1, 0 ?
R2 point(-1, 0);
int nv = mesh.FindVertexNumber(point);

Related topics :

SetGeometryOrder

Location :

MeshBase.hxx
MeshBase.cxx

GetReferenceElement

Syntax

ElementGeomReference<Dimension2> GetReferenceElement(int i) const
Vector<ElementGeomReference<Dimension2>* > GetReferenceElement() const

This method returns the finite element used for the i-th element. In the second syntax, you can also retrieve the finite elements for all the elements.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// finite element for element 4
ElementGeomReference<Dimension2> Fb = mesh.GetReferenceElement(4);

Related topics :

SetGeometryOrder

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetAverageLengthOnBoundary

Syntax

Real_wp GetAverageLengthOnBoundary(const Mesh<Dimension2>& mesh) const

This method returns the mesh size for all the edges of the outer boundary of the PML region.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// mesh size on the outer boundary of PML
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0);

double h = pml.GetAverageLengthOnBoundary(mesh);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

IsPointOnBoundaryRef

Syntax

bool IsPointOnBoundaryRef(int ne, R2 point) const

This method returns true if point is located on the referenced edge ne.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// point on edge 13 ?
int ne = 13;
R2 point(0.5, -2.6);
if (mesh.IsPointOnBoundaryRef(ne, point))
  cout << "Point on edge " << endl;

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

GetAllReferencesAroundPoint

Syntax

void GetAllReferencesAroundPoint(int nv, IVect& ref) const

This method returns all the references of the edges sharing the vertex nv.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// references of neighboring edges ?
int nv = 24;
IVect ref;
mesh.GetAllReferencesAroundPoint(nv, ref);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

GetPointWithOtherReference

Syntax

int GetPointWithOtherReference(int nv, int ref, R2& other_pt, int& other_ref) const

This method returns the number of a vertex close to vertex nv. These two vertices are on a same edge whose reference is different from ref. other_pt are the coordinates of this other vertex, and other_ref the other reference.

Parameters

nv (in)
vertex number
ref (in)
reference associated with vertex nv
other_pt (out)
coordinates of the other point found
other_ref (out)
other reference found

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// around vertex 24, we have already an edge of reference 2
int nv = 24; int ref = 2;
// and we want to get another point belonging to an edge with a
// different reference (different from 2)
int other_ref; R2 other_pt;
int nv2 = mesh.GetPointWithOtherReference(nv, ref, other_pt, other_ref);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

ClearCurves

Syntax

void ClearCurves()

This method clears informations about predefined curves.

Example :

Mesh<Dimension2> mesh;

// adding a predefined curve
mesh.SetCurveType(1, mesh.CURVE_ELLIPSE);

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// clearing predefined curves
mesh.ClearCurves();

// reading another mesh (without ellipse)
mesh.Read("test2.mesh");

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundariesInline.cxx

RemoveReference

Syntax

void RemoveReference(int ref)

This method removes all the edges whose reference is ref.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// removing edges of reference 3
mesh.RemoveReference(3);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

SortPeriodicBoundaries

Syntax

void SortPeriodicBoundaries()

This method changes vertex numbers of the mesh so that matching edges (of two opposite boundaries) have the same orientation. This method is automatically called when reading a mesh, so you should not call this method.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// if you set a periodic condition after reading
TinyVector<int, 2> num;
num(0) = 3; num(1) = 5;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// changing vertex numbers because of that periodic condition
mesh.SortPeriodicBoundaries();
// then you need to call those two functions
mesh.SortBoundariesRef();
mesh.FindConnectivity();

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

SortBoundariesRef

Syntax

void SortBoundariesRef()

This method sorts referenced edges so that the first extremity has the lowest number. This method is actually called when reading a mesh, so you don't need to call that method in regular use.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// if you set a periodic condition after reading
TinyVector<int, 2> num;
num(0) = 3; num(1) = 5;
mesh.AddPeriodicCondition(num, BoundaryConditionEnum::PERIODIC_CTE);

// changing vertex numbers because of that periodic condition
mesh.SortPeriodicBoundaries();
// then you need to call those two functions
mesh.SortBoundariesRef();
mesh.FindConnectivity();

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

ConstructCrossReferenceBoundaryRef

Syntax

void ConstructCrossReferenceBoundaryRef()

This method links referenced edges with global edges. This should not be called in regular use.

Example :

Mesh<Dimension2> mesh;

mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// already done, so useless to call it
mesh.ConstructCrossReferenceBoundaryRef();

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

AddBoundaryEdges

Syntax

void AddBoundaryEdges()

This method detects isolated edges (located on the boundary) without reference, and gives them a new reference. It is already called when reading a mesh. It can be useful to call this function if you constructed manually a mesh.

Example :

Mesh<Dimension2> mesh;

// you construct manually a mesh
mesh.ReallocateVertices(5);
mesh.Vertex(0).Init(1.0, 2.4);
// ...
mesh. ReallocateElements(3);
mesh.Element(0).InitTriangular(0, 2, 5, 1);

// you call FindConnectivity to construct all the edges
mesh.FindConnectivity();

// you did not provide boundary edges in the array BoundaryRef
// a way to retrieve them is to call AddBoundaryEdges
mesh.AddBoundaryEdges();

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

ProjectPointsOnCurves

Syntax

void ProjectPointsOnCurves()

This method projects intermediary points of edges on predefined curves of the mesh. In regular use, you don't have to call that method.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// SetCurve should have been called before but you call it after 
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);

// in that case, you need to call ProjectPointsOnCurves to project points on the curved surfaces
mesh.ProjectPointsOnCurves();

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

RedistributeReferences

Syntax

void RedistributeReferences()

This method redistribute references for boundary edges. For each isolated curve, a different reference is given. This method is useful if you have a mesh without referenced edges, since each boundary will be attributed with a different reference.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// the mesh does not have referenced edges, you call RedistributeReferences
// to distinguish the different boundaries/holes of the mesh
mesh.RedistributeReferences();

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

CheckCurveParameter

Syntax

void CheckCurveParameter(int ref)

This method checks that the parameters associated with the reference ref are correct. If one parameter is not correct, a message and displayed and the program is aborted.

Example :

Mesh<Dimension2> mesh;

// setting curve parameters
mesh.SetCurveType(1, mesh.CURVE_CIRCLE);
// parameters : center of circle and radius
Vector<double> param(3);
param(0) = 0;
param(1) = 0;
param(2) = 2.0;
mesh.SetCurveParameter(1, param);

// you can check that parameters of reference 1 are correct
mesh.CheckCurveParameter(1);

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

ReadCurveType

Syntax

void ReadCurveType(const IVect& ref, string keyword)

This method sets curve types of all the references in array ref. The curve type is provided through the string keyword. The method SetCurveType should be preferred.

Example :

Mesh<Dimension2> mesh;

// predefined curves
IVect ref(3);
ref(0) = 1; ref(1) = 4; ref(2) = 5;
mesh.ReadCurveType(ref, string("CIRCLE"));

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

ProjectToCurve

Syntax

void ProjectToCurve(R2& ptM, int ref) const

This method projects a point to the curve whose reference is ref. If this reference is not associated with a predefined curve, the point will not be modified.

Example :

Mesh<Dimension2> mesh;

// predefined curves
mesh.SetCurveType(2, mesh.CURVE_CIRCLE);

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// then you can project a point to the curve 2 (which is a circle)
R2 pt(0.3, 0.5);
mesh.ProjectToCurve(pt, 2);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetPointsOnCurve

Syntax

void GetPointsOnCurve(int i, int ref, const VectReal_wp& s, VectR2& pts) const

This method computes intermediary points of a given edge.

Parameters

i (in)
edge number
ref (in)
reference of the edge
s (in)
local coordinates on interval [0, 1]
pts (out)
computed points on the edge

Example :

Mesh<Dimension2> mesh;

// predefined curves
mesh.SetCurveType(2, mesh.CURVE_CIRCLE);

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// then you can compute intermediary points of an edge
int num_edge = 10;
int ref = mesh.BoundaryRef(num_edge).GetReference();
// you want intermediary points for which values of parameter s ?
// values should be on interval [0, 1]
// s = 0 -> first extremity of edge, s = 1 -> second extremity
VectReal_wp s(3); VectR2 points(3);
s(0) = 0.14; s(1) = 0.46; s(2) = 0.82;
mesh.GetPointsOnCurve(num_edge, ref, s, points);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetAllPointsOnReference

Syntax

void GetAllPointsOnReference(int ref, VectR2& points, VectR2& normale)
void GetAllPointsOnReference(int ref_target, IVect ref_cond, VectR2& points, VectR2& normale)

This method retrieves all the points and normales for a curve whose reference is ref. In the second syntax, we consider all the references i such that ref_cond(i) is equal to ref_target.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(3);
mesh.Read("test.mesh");

// points and normale of reference 3
VectR2 points, normale;
mesh.GetAllPointsOnReference(3, points, normale);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetDistantPointsOnReference

Syntax

void GetDistantPointsOnReference(int ref, VectR2& points, VectR2& normale, int nb_points)
void GetDistantPointsOnReference(int ref_target, IVect ref_cond, VectR2& points, VectR2& normale, int nb_points)

This method retrieves some points and normales for a curve whose reference is ref. In the second syntax, we consider all the references i such that ref_cond(i) is equal to ref_target. It tries to find nb_points differents points which are distant enough (possibly non-collinear).

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// 3 distants points and normale of reference 5
VectR2 points, normale;
mesh.GetDistantPointsOnReference(5, points, normale, 3);

Related topics :

SetGeometryOrder

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersPlane

Syntax

void FindParametersPlane(int ref, VectReal_wp& param)

This method finds coefficients a, b, c (equation ax + by + c = 0) of a line if we assume that edges with reference equal to ref are on this line.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// coefficients a, b, c of line whose reference is 4 
VectReal_wp param(3);
mesh.FindParametersPlane(4, param);

Related topics :

FindParametersCurve

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersSpline

Syntax

void FindParametersSpline(int ref, VectReal_wp& param)

This method finds parameters of a cubic spline if we assume that edges with reference equal to ref are on this spline.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// coefficients of a spline
VectReal_wp param;
mesh.FindParametersSpline(4, param);

Related topics :

FindParametersCurve

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersCircle

Syntax

void FindParametersCircle(int ref, VectReal_wp& param)

This method finds parameters xc, yc, r (equation (x-xc)2 + (y-yc)2 = r2) of a circle if we assume that edges with reference equal to ref are on this circle.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// parameters of a circle whose reference is 4 
VectReal_wp param(3);
mesh.FindParametersCircle(4, param);

Related topics :

FindParametersCurve

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersEllipse

Syntax

void FindParametersEllipse(int ref, VectReal_wp& param)

This method finds parameters xc, yc, a, b of an ellipse if we assume that edges with reference equal to ref are on this ellipse. The considered equation is

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// parameters of an ellipse whose reference is 4 
VectReal_wp param(4);
mesh.FindParametersEllipse(4, param);

Related topics :

FindParametersCurve

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersPeanut

Syntax

void FindParametersPeanut(int ref, VectReal_wp& param)

This method finds parameters xc, yc, a, b, c of a peanut if we assume that edges with reference equal to ref are on this peanut. The considered equation is

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// parameters of an ellipse whose reference is 4 
VectReal_wp param(4);
mesh.FindParametersPeanut(4, param);

Related topics :

FindParametersCurve

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FindParametersCurve

Syntax

void FindParametersCurve()

This method finds parameters of all predefined curves of a mesh. It is automatically called when reading a mesh, so you don't need to call it.

Example :

Mesh<Dimension2> mesh;

// reading mesh
mesh.SetGeometryOrder(4);
mesh.Read("test.mesh");

// useless (already called in Read)
mesh.FindParametersCurve();

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetCenterRadialPML

Syntax

R2 GetCenterRadialPML()

This method returns the center of radial PMLs.

Example :

Mesh<Dimension2> mesh;

R2 center = mesh.GetCenterRadialPML();

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

AddPMLElements

Syntax

void AddPMLElements(int nb_layers, int num, const Mesh<Dimension2>& mesh)

This method adds PML layers to a mesh. You need to provide the number of layers. If you want to compute automatically this number (from the mesh size and thickness of PML), the method AddPML should be used.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0); // number of the region is 0

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(mesh.PML_NEGATIVE_SIDE, mesh.PML_BOTH_SIDES, mesh.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// mesh is modified in order to add those PML
pml.AddPMLElements(nb_layers, 0, mesh);

Related topics :

SetAdditionPML

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetNbParameters

Syntax

int GetNbParameters() const

This method returns the number of parameters associated with the PML region.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0); // number of the region is 0

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(mesh.PML_NEGATIVE_SIDE, mesh.PML_BOTH_SIDES, mesh.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// number of parameters ?
int nb_param = pml.GetNbParameters();

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

FillParameters

Syntax

void FillParameters(VectReal_wp& all_param, int& param_num)

This method appends to the array all_param the parameters associated with the PML region. The values are appended from the index param_num which is incremented.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0); // number of the region is 0

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(mesh.PML_NEGATIVE_SIDE, mesh.PML_BOTH_SIDES, mesh.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// storing parameters of the PML
VectReal_wp all_param(20);
int num_param = 10;
pml.FillParameters(all_param, num_param);
// all_param(10:10+num_param) will contain parameters of the PML

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

SetRegion

Syntax

void SetRegion(const VectReal_wp& all_param, int& param_num)

This method sets parameters of the PML region with the array all_param. The values are read from the index param_num which is incremented. This method is the reciprocal of FillParameters.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// if you want only one PML region
mesh.ReallocatePmlAreas(1);
PmlRegionParameter<Dimension2>& pml = mesh.GetPmlArea(0); // number of the region is 0

// adding PML for x < xmin, y < ymin and y > ymax
int nb_layers = 2;
pml.SetAdditionPML(mesh.PML_NEGATIVE_SIDE, mesh.PML_BOTH_SIDES, mesh.PML_NO, nb_layers);
pml.SetThicknessPML(1.5);

// storing parameters of the PML
VectReal_wp all_param(20);
int num_param = 10;
pml.FillParameters(all_param, num_param);
// all_param(10:10+num_param) will contain parameters of the PML

// you can set a region with the stored parameters
PmlRegionParameter<Dimension2> pml2;
num_param = 10;
pml2.SetRegion(all_param, num_param);

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GetBoundaryMesh

Syntax

void GetBoundaryMesh(int ref, SurfacicMesh<Dimension2>& mesh, const IVect& ref_cond)
void GetBoundaryMesh(int ref, Mesh<Dimension2>& mesh, IVect& ListeEdge, IVect& ListeVertex, IVect& NumElem, IVect& NumLoc, const IVect& ref_cond)

This method extracts all the referenced edges whose reference i satisfies ref_cond(i) = ref.

Parameters

ref (in)
reference target
mesh (out)
boundary mesh created
ref_cond (in)
we extract all edges of reference i such that ref_cond(i) = ref
ListeEdge (out)
extracted edge numbers
ListeVertex (out)
extracted vertex numbers
NumElem (out)
numbers of elements that own the extracted edges
NumLoc (out)
local positions of the extracted edges in their elements

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// extracting reference 2 and 3
IVect ref_cond(mesh.GetNbReferences()+1);
ref_cond.Fill(0);
ref_cond(2) = 1; ref_cond(3) = 1;
SurfacicMesh<Dimension2> surf_mesh;
mesh.GetBoundaryMesh(1, surf_mesh, ref_cond);

// you can also retrieve a simple 2-D mesh
Mesh<Dimension2> sub_mesh;
// and you get additional arrays giving vertex numbers, edge numbers,
// and for each extracted edge, the element and local position on the original mesh
IVect ListeEdge, ListeVertex, NumElem, NumLocalBoundary;
mesh.GetBoundaryMesh(1, sub_mesh, ListeEdge, ListeVertex, NumElem,
                     NumLocalBoundary, ref_cond);

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

GenerateSurfaceOfRevolution

Syntax

void GenerateSurfaceOfRevolution(int ref, Mesh<Dimension3>& mesh, IVect& FaceToEdge, IVect& Vertex, VectReal_wp& AngleVertex, const IVect& ref_cond)

This method generates a triangular mesh of an axisymmetric surface by rotation of a 2-D curve around axis Oz, as shown in figure below.

Parameters

ref (in)
reference target
mesh (out)
3-D created mesh by extrusion
ref_cond (in)
we extrude all edges of reference i such that ref_cond(i) = ref
FaceToEdge (out)
for each 3-D face i, FaceToEdge(i) is the number of the extruded edge
Vertex (out)
for each 3-D vertex i, Vertex(i) is the number of the extruded 2-D vertex
AngleVertex (out)
AngleVertex(i) is the angle θ for the vertex i
Loading     ⇒     Loading

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

// 2-D curve consists of edges with reference 2 and 3
IVect ref_cond(mesh.GetNbReferences()+1);
ref_cond.Fill(0);
ref_cond(2) = 1;
ref_cond(3) = 1;

// generation of 3-D triangular mesh
Mesh<Dimension2> mesh3d;
IVect FaceToEdge, VertexToVert; VectReal_wp AngleVertex;
mesh.GenerateSurfaceOfRevolution(1, mesh3d, FaceToEdge, VertexToVert, AngleVertex, ref_cond);

// then you can write the result mesh
mesh3d.Write("axi.mesh");

Related topics :

SetAdditionPML

Location :

Mesh2DBoundaries.hxx
Mesh2DBoundaries.cxx

TreatPeriodicCondition

Syntax

void TreatPeriodicCondition()

This method finds periodic dofs, it is called in method NumberMesh, so you should not call that method in regular use.

Example :

Mesh<Dimension2> mesh;

// reading the mesh
mesh.Read("exemple.mesh");

MeshNumbering<Dimension2> mesh_num(mesh);

// if you want to treat periodic condition apart
mesh_num.treat_periodic_condition_during_number = false;
mesh_num.NumberMesh();

// then you compute the periodic dofs
mesh_num.TreatPeriodicCondition();

Related topics :

AddPeriodicCondition

Location :

NumberMesh.hxx
NumberMesh.cxx

GetNbElt

Syntax

int GetNbElt() const

This method returns the number of elements in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we display the number of elements in the mesh
cout << "Number of elements " << mesh.GetNbElt();

Related topics :

Element

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetNbVertices

Syntax

int GetNbVertices() const

This method returns the number of vertices in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we display the number of vertices in the mesh
cout << "Number of vertices " << mesh.GetNbVertices();

Related topics :

Vertex

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetNbEdges

Syntax

int GetNbEdges() const

This method returns the number of edges in the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we display the number of edges in the mesh
cout << "Number of edges " << mesh.GetNbEdges();

Related topics :

GetEdge

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetOriginalNeighborReference

Syntax

const IVect& GetOriginalNeighborReference() const

This method returns the list of original references of the mesh. This method is relevent only in parallel computation when the mesh is split. Because of this splitting, new references are created (for the interface between two processors) that can coincide with previous references. So you can know the original reference number for a given referenced edge.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");

// we assume the mesh is split

// then you can access the original references
const IVect& original_ref = mesh.GetOriginalNeighborReference();

// original reference of reference 8 :
int ref0 = original_ref(8);

Location :

MeshBase.hxx
MeshBaseInline.cxx

SetOriginalNeighborReference

Syntax

void SetOriginalNeighborReference(const IVect& ref)

This method sets the list of original references of the mesh. This method is relevent only in parallel computation when the mesh is split. Because of this splitting, new references are created (for the interface between two processors) that can coincide with previous references. This method is called during the splitting process.

Related topics :

GetOriginalNeighborReference

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetPathName

Syntax

string GetPathName() const

This method returns path of the directory where meshes are searched when ConstructMesh is called.

Example :

Mesh<Dimension2> mesh;
// path ?
cout << "Path where meshes are searched " << mesh.GetPathName() << endl;

Related topics :

SetPathName

Location :

MeshBase.hxx
MeshBaseInline.cxx

SetPathName

Syntax

void SetPathName(const string& path)

This method sets the path of the directory where meshes are searched when ConstructMesh is called.

Example :

Mesh<Dimension2> mesh;
// path ?
mesh.SetPathName("MAILLAGES/");

Related topics :

GetPathName

Location :

MeshBase.hxx
MeshBaseInline.cxx

IsElementAffine

Syntax

bool IsElementAffine(int i) const

This method returns true if the element i is affine. An affine element is an element such that the transformation Fi (that transforms the reference element to an element of the mesh) is affine. For example, a parallelogram is affine whereas a trapezoid is not affine. Straight triangles and tetrahedra are affine.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// then you can see if an element is affine
for (int i = 0; i < mesh.GetNbElt(); i++)
  if (mesh.IsElementAffine(i))
    cout << "Element " << i << " affine ";

Related topics :

ReorientElements

Location :

MeshBase.hxx
MeshBaseInline.cxx

ReallocateElements

Syntax

void ReallocateElements(int n)

This method changes the number of elements present in the mesh. Previous elements stored will be lost, so you have to construct all the elements (with InitTriangular, InitQuadrangular for instance).

Example :

Mesh<Dimension2> mesh;
// we construct a mesh with two elements
mesh.ReallocateElements(2);
int ref = 1;
mesh.Element(0).InitTriangular(0, 1, 5, ref);
mesh.Element(1).InitTriangular(1, 2, 5, ref);

Related topics :

Element

Location :

MeshBase.hxx
MeshBaseInline.cxx

ResizeElements

Syntax

void ResizeElements(int n)

This method changes the number of elements present in the mesh. Previous elements stored will be conserved, so that you have to construct only the new elements (if there are more elements than before).

Example :

Mesh<Dimension2> mesh;
// we construct a mesh with two elements
mesh.ReallocateElements(2);
int ref = 1;
mesh.Element(0).InitTriangular(0, 1, 5, ref);
mesh.Element(1).InitTriangular(1, 2, 5, ref);
// then we add two new elements
mesh.ResizeElements(4);
IVect num(3); num(0) = 2; num(1) = 4; num(2) = 5;
mesh.Element(2).Init(num, ref);
num(0) = 2; num(1) = 4; num(2) = 3;
mesh.Element(3).Init(num, ref);

Related topics :

Element

Location :

MeshBase.hxx
MeshBaseInline.cxx

ClearElements

Syntax

void ClearElements()

This method removes all the elements present in the mesh

Example :

Mesh<Dimension2> mesh;
// we construct a mesh with two elements
mesh.ReallocateElements(2);
int ref = 1;
mesh.Element(0).InitTriangular(0, 1, 5, ref);
mesh.Element(1).InitTriangular(1, 2, 5, ref);
mesh.ResizeElements(4);
IVect num(3); num(0) = 2; num(1) = 4; num(2) = 5;
mesh.Element(2).Init(num, ref);
num(0) = 2; num(1) = 4; num(2) = 3;
mesh.Element(3).Init(num, ref);
// then you clear all those modifications
mesh.ClearElements();

Related topics :

Element

Location :

MeshBase.hxx
MeshBaseInline.cxx

Element

Syntax

Face<Dimension2>& Element(int i)

This method gives access to an element of the mesh.

Example :

Mesh<Dimension2> mesh;
// we construct a mesh with two elements
mesh.ReallocateElements(2);
int ref = 1;
mesh.Element(0).InitTriangular(0, 1, 5, ref);
mesh.Element(1).InitTriangular(1, 2, 5, ref);
// you can also print an element
cout << "Second element : " << mesh.Element(1);

// retrieve the vertex numbers of an element :
int n0 = mesh.Element(1).numVertex(0);
int n1 = mesh.Element(1).numVertex(1);
// etc

Related topics :

ReallocateElements

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetTypeElement

Syntax

int GetTypeElement(int i) const

This method returns the type of the element, i.e. 0 if it is a triangle, 1 for a quadrilateral.

Example :

Mesh<Dimension2> mesh;
// we construct a mesh with two elements
mesh.ReallocateElements(2);
int ref = 1;
mesh.Element(0).InitTriangular(0, 1, 5, ref);
mesh.Element(1).InitQuadrangular(1, 2, 5, 3, ref);
// you can also print the type of the element
cout << "Type of second element : " << mesh.GetTypeElement(1);

Related topics :

Element

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetEdge

Syntax

Edge<Dimension2>& GetEdge(int i)

This method gives access to an edge of the mesh.

Example :

Mesh<Dimension2> mesh;
// we read the mesh
mesh.Read("test.mesh");
// we count edges on the boundary
int nb_edges_bound = 0;
for (int i = 0; i < mesh.GetNbEdges(); i++)
  if (mesh.GetEdge(i).GetNbElements() == 1)
    nb_edges_bound++;

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBaseInline.cxx

ReallocateEdges

Syntax

void ReallocateEdges(int n)

This method changes the number of edges stored in the mesh. Previous edges will be removed. Usually, we don't enter directly edge connectivity, but we call FindConnectivity.

Example :

Mesh<Dimension2> mesh;
// you specify the number of edges
mesh.ReallocateEdges(5);
// then you can set edges
mesh.GetEdge(0).Init(0, 2, 1);
mesh.GetEdge(1).Init(1, 12, 1);
// ...

Related topics :

ResizeEdges

Location :

MeshBase.hxx
MeshBase.cxx

ResizeEdges

Syntax

void ResizeEdges(int n)

This method changes the number of edges stored in the mesh. Previous edges will be kept. Usually, we don't enter directly edge connectivity, but we call FindConnectivity.

Example :

Mesh<Dimension2> mesh;
// you specify the number of edges
mesh.ResizeEdges(5);
// then you can set edges
mesh.GetEdge(0).Init(0, 2, 1);
mesh.GetEdge(1).Init(1, 12, 1);
// ...

Related topics :

ReallocateEdges

Location :

MeshBase.hxx
MeshBase.cxx

ReallocateVertices

Syntax

void ReallocateVertices(int n)

This method changes the number of vertices stored in the mesh. Previous vertices will be removed and replaced by 0. After, you can modify each vertex with the method Vertex.

Example :

Mesh<Dimension2> mesh;
// you specify the number of vertices
mesh.ReallocateVertices(5);
// then you fill each point
mesh.Vertex(0).Init(-1.5, -1.5); mesh.Vertex(1) = R2(0, -1.5);
mesh.Vertex(2).Init(1.5, 0); mesh.Vertex(3) = R2(1.5, 1.5);
mesh.Vertex(4) = R2(0, 1.5);

Related topics :

ResizeVertices

Location :

MeshBase.hxx
MeshBase.cxx

ResizeVertices

Syntax

void ResizeVertices(int n)

This method changes the number of vertices stored in the mesh. Previous vertices will be conserved. After, you have to specify each new vertex with the method Vertex.

Example :

Mesh<Dimension2> mesh;
// you specify the number of vertices
mesh.ReallocateVertices(5);
// then you fill each point
mesh.Vertex(0).Init(-1.5, -1.5); mesh.Vertex(1) = R2(0, -1.5);
mesh.Vertex(2).Init(1.5, 0); mesh.Vertex(3) = R2(1.5, 1.5);
mesh.Vertex(4) = R2(0, 1.5);
// we add a new point
mesh.ResizeVertices(6);
mesh.Vertex(5).Init(0.2, 0.2);

Related topics :

ReallocateVertices

Location :

MeshBase.hxx
MeshBase.cxx

ClearVertices

Syntax

void ClearVertices(int n)

This method removes all the vertices of the mesh.

Example :

Mesh<Dimension2> mesh;
// you specify the number of vertices
mesh.ReallocateVertices(5);
// then you fill each point
mesh.Vertex(0).Init(-1.5, -1.5); mesh.Vertex(1) = R2(0, -1.5);
mesh.Vertex(2).Init(1.5, 0); mesh.Vertex(3) = R2(1.5, 1.5);
mesh.Vertex(4) = R2(0, 1.5);
// you can clear all the vertices
mesh.ClearVertices();

Related topics :

ReallocateVertices

Location :

MeshBase.hxx
MeshBaseInline.cxx

Vertex

Syntax

R2& Vertex(int i)

This method gives access to a vertex of the mesh.

Example :

Mesh<Dimension2> mesh;
// you specify the number of vertices
mesh.ReallocateVertices(5);
// then you fill each point
mesh.Vertex(0).Init(-1.5, -1.5); mesh.Vertex(1) = R2(0, -1.5);
mesh.Vertex(2).Init(1.5, 0); mesh.Vertex(3) = R2(1.5, 1.5);
mesh.Vertex(4) = R2(0, 1.5);

Related topics :

ReallocateVertices

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetXmin, GetXmax, GetYmin, GetYmax

Syntax

Real_wp GetXmin()
Real_wp GetYmin()
Real_wp GetXmax()
Real_wp GetYmax()

These methods return the minimum/maximum of x or y-coordinate of vertices in the mesh. These values are computed when you read or construct the mesh.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");
cout << "Mesh is contained in box [" << mesh.GeXmin() << ", " << mesh.GetXmax() << "] x [ "
<< mesh.GetYmin() << ", " << mesh.GetYmax() << "] ";

// if you construct the mesh manually
mesh.Clear();
// you specify the number of vertices
mesh.ReallocateVertices(5);
// then you fill each point
mesh.Vertex(0).Init(-1.5, -1.5); mesh.Vertex(1) = R2(0, -1.5);
mesh.Vertex(2).Init(1.5, 0); mesh.Vertex(3) = R2(1.5, 1.5);
mesh.Vertex(4) = R2(0, 1.5);
// the bounding box is updated when you call ReorientElements
mesh.ReorientElements();
cout << "Mesh is contained in box [" << mesh.GeXmin() << ", " << mesh.GetXmax() << "] x [ "
<< mesh.GetYmin() << ", " << mesh.GetYmax() << "] ";

Location :

MeshBase.hxx
MeshBaseInline.cxx

GetVerticesElement

Syntax

void GetVerticesElement(int i, VectR2& vertices)

This method retrieves the vertices of an element.

Example :

Mesh<Dimension2> mesh;
// mesh is read on a file
mesh.Read("test.mesh");
VectR2 s;
// loop over elements
for (int i = 0; i < mesh.GetNbElt(); i++)
  {
    mesh.GetVerticesElement(i, s);
    // then you can perform operations over vertices s
    cout << "Second vertex " << s(1) << endl;
  }

Related topics :

Vertex

Location :

MeshBase.hxx
MeshBase.cxx

GetDofsElement

Syntax

void GetDofsElement(int i, VectR2& points, const ElementGeomReference<Dimension2>& fem)

This method retrieves the locations of degrees of freedom on an element. It does not need a numbering of the mesh, but a finite element object fem where points on degrees of freedom are defined. The finite element must be of the same type as the element (e.g. triangular if the element is a triangle).

Example :

Mesh<Dimension2> mesh;

// mesh is read on a file
mesh.Read("test.mesh");

// finite element with dofs
QuadrangleLobatto quad;
quad.ConstructFiniteElement(5);

// Gauss-Lobatto points after transformation Fi for a given element i
VectR2 points; int i = 10;
mesh.GetDofsElement(i, points, quad);
DISP(points);

Location :

MeshBase.hxx
MeshBase.cxx

SetInputData

Syntax

void SetInputData(const string& keyword, VectString& param)

This method is used during the reading of the data file. You can look at the section devoted to the data file and see the fields related to the mesh.

Example :

Mesh<Dimension2> mesh;
// some parameters are read on a data file
ReadInputFile("test.ini", &mesh);

Related topics :

Data files

Location :

MeshBase.hxx
MeshBase.cxx

PointsOnSameEdge

Syntax

bool PointsOnSameEdge(int n1, int n2, int& num_edge)

This method returns true if vertices n1 and n2 are the extremities of the same edge. On exit, num_edge contains the edge number.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// we test if two points share the same edge
int n1 = 2, n2 = 5, ne = -1;
if (PointsOnSameEdge(n1, n2, ne))
  cout << "Points are on the edge " << ne << endl;

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

GetNeighboringElementsAroundVertex

Syntax

void GetNeighboringElementsAroundVertex(int i, IVect& num)

This method fills the array num with the numbers of the elements sharing the point i.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// we find all the elements around vertex 3
IVect num;
mesh.GetNeighboringElementsAroundVertex(3, num);
cout << "Elements around vertex 3 : " << num;

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

GetNeighboringEdgesAroundVertex

Syntax

void GetNeighboringEdgesAroundVertex(int i, IVect& num)

This method fills the array num with the numbers of the edges sharing the point i.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// we find all the edges around vertex 3
IVect num;
mesh.GetNeighboringEdgesAroundVertex(3, num);
cout << "Edges around vertex 3 : " << num;

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

SetVerticesToBeAdded

Syntax

void SetVerticesToBeAdded(const VectR2& pts)

This method sets the points that will be added when ConstructMesh will be called.

Example :

Mesh<Dimension2> mesh;

// you set the vertices to add to the final mesh
VectR2  s(2);
s(0).Init(-2.0, -0.8); s(1).Init(-1.0, -0.5);
mesh.SetVerticesToBeAdded(s);

// then ConstructMesh will add these vertices
VectString param(1);
param(0) = string("carre.mesh");
mesh.ConstructMesh(mesh.HYBRID_MESH, param);

Related topics :

ConstructMesh

Location :

MeshBase.hxx
MeshBase.cxx

SetVerticesToBeRefined

Syntax

void SetVerticesToBeRefined(const VectR2& pts, int lvl, Real_wp ratio)

This method sets the points that will be refined when ConstructMesh will be called.

Example :

Mesh<Dimension2> mesh;

// you set the vertices around which the final mesh will be refined
// level : refinement level, ratio : refinement ratio
int level = 2; Real_wp ratio = 3.0;
VectR2  s(2);
s(0).Init(-2.0, -0.8); s(1).Init(-1.0, -0.5);
mesh.SetVerticesToBeRefined(s, level, ratio);

// then ConstructMesh will refine locally the mesh around these vertices
VectString param(1);
param(0) = string("carre.mesh");
mesh.ConstructMesh(mesh.HYBRID_MESH, param);

Related topics :

ConstructMesh

Location :

MeshBase.hxx
MeshBase.cxx

CreateReferenceVertices

Syntax

void CreateReferenceVertices(IVect& ref_vertex) const

This method attributes a reference for each vertex. The reference is 0 if the vertex does not belong to a referenced edge, otherwise the reference is equal to the reference of the edge.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
IVect ref_vertex;
mesh.CreateReferenceVertices(ref_vertex);

Related topics :

BoundaryRef

Location :

MeshBase.hxx
MeshBase.cxx

GetMeshSize

Syntax

Real_wp GetMeshSize() const

This method returns the mesh size, which is evaluated as the average value of edge length.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// we return the mesh size
cout << "h = " << mesh.GetMeshSize();

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

ClearConnectivity

Syntax

void ClearConnectivity()

This method clears edge numbers for each face, and all the edges. Therefore, the mesh is stored as if you didn't call the method FindConnectivity (which is called when reading a mesh, and lots of other methods). So, it can reduce the storage required by the mesh. Some methods may disfunction when clearing connectivity.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// then you clear arrays containing edges
mesh.ClearConnectivity();

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

Clear

Syntax

void Clear()

This method clears the mesh.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// then you can clear all the mesh
mesh.Clear();

Related topics :

FindConnectivity

Location :

MeshBase.hxx
MeshBase.cxx

ClearInput

Syntax

void ClearInput()

This method clears the input parameters of the mesh. These parameters can have been set by calling SetInputData for example.

Example :

Mesh<Dimension2> mesh;
// if you called SetInputData

// you can clear input data
mesh.ClearInput();

Location :

MeshBase.hxx
MeshBase.cxx

AddOverlapRegion

Syntax

void AddOverlapRegion(int n, int& nb_vert, int& nb_elt, Vector<bool>& VertexOnSubdomain, Vector<bool>& ElementOnSubdomain) const

This method adds n layers around a subdomain of the mesh. It is used by the method SplitIntoBoxes, SplitMetis and equivalents if overlapping subdomains are required by the user.

Parameters

n(in)
number of layers to add
nb_vert(in/out)
number of vertices in the subdomain (this number is updated)
nb_elt(in/out)
number of elements in the subdomain (this number if updated)
VertexOnSubdomain(in/out)
if VertexOnSubdomain(i) is true the vertex i lies in the subdomain
ElementOnSubdomain(in/out)
if ElementOnSubdomain(i) is true the element i lies in the subdomain

Example :

Mesh<Dimension2> mesh;

// we read mesh
mesh.Read("test.mesh");

// for example you consider a single element
int nb_elt = 1; int num_elem = 23;
Vector<bool> ElementOnSubdomain(mesh.GetNbElt()); ElementOnSubdomain.Fill(false);
ElementOnSubdomain(num_elem) = true;
// and vertices of this element
int nb_vert = mesh.Element(num_elem).GetNbVertices();
Vector<bool> VertexOnSubdomain(mesh.GetNbVertices()); VertexOnSubdomain.Fill(false);
for (int i = 0; i < nb_vert; i++)
  VertexOnSubdomain(mesh.Element(num_elem).numVertex(i)) = true;

// then you add 3 layers of elements around this initial element
mesh.AddOverlapRegion(3, nb_vert, nb_elt, VertexOnSubdomain, ElementOnSubdomain);

// you can create a mesh with the updated informations
mesh.CreateSubMesh(sub_mesh, nb_vert, nb_elt, VertexOnSubdomain, ElementOnSubdomain);
sub_mesh.Write("subdomain.mesh");

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

SplitIntoBoxes

Syntax

void SplitIntoBoxes(int nx, int ny, int nz, int& nb_parts, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh within a cartesian grid, with nx cells along x, ny cells along y. The number of subdomains is returned in variable nb_parts, and can be different from nx*ny*nz, if there are holes for example. You can define an overlapping level by increasing noverlap, if equal to 0, there will be no overlapping. On exit, num_elem contains the element numbers for each subdomain.

Parameters

nx(in)
number of subdivisions along x
ny(in)
number of subdivisions along y
nz(in)
number of subdivisions along z
nb_parts(in/out)
number of subdomains
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");
// we split mesh into 2x2 grid without overlap
int nb_parts = 0;
Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
mesh.SplitIntoBoxes(2, 2, 1, nb_parts, sub_mesh, 0);
// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

SplitMetis

Syntax

void SplitMetis(int nb_parts, IVect& weight_elt, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh into nb_parts balanced parts by using Metis routines. weight_elt denotes the weight attributed to each element. You can define an overlapping level by increasing noverlap, if equal to 0, there will be no overlapping. On exit, num_elem contains the element numbers for each subdomain.

Parameters

nb_parts(in)
number of subdomains
weight_elt(in)
weight for each element
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// we split mesh into 7 parts with weight equal to 1
IVect weight_elt(mesh.GetNbElt()); weight_elt.Fill(1);
Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
mesh.SplitMetis(7, weight_elt, num_elt, sub_mesh, 0);

// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitIntoBoxes

Location :

MeshBase.hxx
MeshBase.cxx

SplitScotch

Syntax

void SplitScotch(int nb_parts, IVect& weight_elt, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh into nb_parts balanced parts by using Scotch routines. weight_elt denotes the weight attributed to each element. You can define an overlapping level by increasing noverlap, if equal to 0, there will be no overlapping. On exit, num_elem contains the element numbers for each subdomain.

Parameters

nb_parts(in)
number of subdomains
weight_elt(in)
weight for each element
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;

// we read mesh
mesh.Read("test.mesh");

// we split mesh into 7 parts with weight equal to 1
IVect weight_elt(mesh.GetNbElt()); weight_elt.Fill(1);
Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
mesh.SplitScotch(7, weight_elt, num_elt, sub_mesh, 0);

// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

SplitConcentric

Syntax

void SplitConcentric(int nb_parts, VectReal_wp& radius, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh by using the radii provided by the user. The subdomain 0 will comprise elements such that r ≤ r0 (here r stands for the radius of polar coordinates), the subdomain 1 will comprise elements such that r0 < r ≤ r1, etc. The last subdomain comprises elements such taht r ≥ rN.

Parameters

nb_parts(in)
number of subdomains
radius(in)
list of radii used to distinguish subdomains
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// we split mesh into 3 concentric parts
int nb_parts = 3;
VectReal_wp radius(2); radius(0) = 0.5; radius(1) = 2.0;

Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
mesh.SplitConcentric(nb_parts, radius, num_elt, sub_mesh, 0);

// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

SplitLayered

Syntax

void SplitLayered(int nb_parts, const Vector<IVect>& ref_layer, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh by using the references provided by the user. The subdomain i will comprise elements whose reference is contained in the array ref_layer(i).

Parameters

nb_parts(in)
number of subdomains
ref_layer(in)
list of references for each subdomain
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// we split mesh into 3 parts
// first part will comprise elements of reference 1,
// second part elements of reference 2 and 3
// last part will contain elements of reference 2
int nb_parts = 3;
Vector<IVect> ref_layer(3);
ref_layer(0).PushBack(1);
ref_layer(1).PushBack(3); ref_layer(1).PushBack(4);
ref_layer(2).PushBack(2);

Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
mesh.SplitLayered(nb_parts, ref_layer, num_elt, sub_mesh, 0);

// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

PartMesh

Syntax

void PartMesh(int nb_parts, const IVect& epart, Vector<IVect>& num_elem, Vector<Mesh<Dimension2> >& sub_mesh, int noverlap) const

This method splits the current mesh by using the array epart provided by the user. The subdomain i will comprise elements j such that epart(j) is equal to i.

Parameters

nb_parts(in)
number of subdomains
epart(in)
subdomain number for each element of the mesh
num_elem(out)
Element numbers for all subdomains
sub_mesh(out)
Meshes of the subdomains
noverlap(in)
overlapping level

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// the user gives the subdomain number in epart
IVect epart; int nb_parts = 4;
epart.Read("epart.dat");

Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
PartMesh(nb_parts, epart, num_elt, sub_mesh, 0);

// then you can write the small meshes
for (int i = 0; i < nb_parts; i++)
  sub_mesh(i).Write(string("subdomain")+to_str(i)+".mesh"); 

Related topics :

SplitMetis

Location :

MeshBase.hxx
MeshBase.cxx

UpdatePeriodicReferenceSplit

Syntax

void UpdatePeriodicReferenceSplit(const IVect& epart, const IVect& NumLoc, Vector<Mesh<Dimension2> >& sub_mesh)

This methods updates references if periodic conditions have been set on the mesh. If the splitting of the mesh has been constructed such that periodic boundaries are separated (one element on a processor, and the opposite element on another processor), a new reference is created for this boundary.

Parameters

epart(in)
subdomain number for each element of the mesh
NumLoc(in)
local number of element i in its subdomain
sub_mesh(in/out)
Meshes of the subdomains

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// the user gives the subdomain number in epart
IVect epart; int nb_parts = 4;
epart.Read("epart.dat");

Vector<IVect> num_elt;
Vector<Mesh<Dimension2> > sub_mesh;
PartMesh(nb_parts, epart, num_elt, sub_mesh, 0);

// numloc(i) is the local element number in the corresponding subdomain
IVect numloc(mesh.GetNbElt());
IVect num(nb_parts); num.Zero();
for (int i = 0; i < mesh.GetNbElt(); i++)
  numloc(i) = num(epart(i))++;

// updating periodic references
UpdatePeriodicReferenceSplit(epart, numloc, sub_mesh);

Location :

MeshBase.hxx
MeshBase.cxx

WriteOrder

Syntax

void WriteOrder(const string& file_name, const MeshNumbering<Dimension2>& mesh_num)

This method writes a mesh in a file, on which references attributed to elements are equal to the order (given by mesh_num). The mesh is not modified when calling this method.

Example :

Mesh<Dimension2> mesh;
// we read mesh
mesh.Read("test.mesh");

// we compute a numbering
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.ComputeVariableOrder(true);
mesh_num.NumberMesh();

// then you can see the result
mesh.WriteOrder("order.mesh", mesh_num);

Related topics :

Write

Location :

MeshBase.hxx
MeshBase.cxx

ExtrudePMLLayer

Syntax

void ExtrudePMLLayer( Mesh<Dimension2>& mesh, R2 axis, Real_wp scal_max, IVect ref_cond, int ref_target, Real_wp delta, int nb_layers, bool circular)

This function adds layers to a given boundary. If the last argument circular is false, the boundary is assumed to be a slanted line (axis is the normale). If the last argument is true, the boundary is assumed to be a circle (axis is the center of the circle). This operation is a typical operation needed to add PMLs to a computational domain). If ref_target is different from -1, the boundary is found by searching edges such that ref_cond(ref) is equal to ref_target.

Example :

// you start from an initial mesh
Mesh<Dimension2> mesh;
mesh.Read("example.mesh");

// thickness of PMLs
Real_wp delta = 0.5;
// number of cells in direction of extrusion
int nb_layers = 2;

// for a slanted boundary, if you want to add layers (such as PMLs)
// you have to give the axis of extrusion (direction where layers will be added)
Real_wp teta = pi_wp/5;
R2 axis(cos(teta), sin(teta));

// you can provide the reference of the boundary where layers will be appended
int my_provided_ref = 2, ref_target = 1;
Vector<int> ref_cond(mesh.GetNbReferences()+1); ref_cond.Zero();
ref_cond(my_provided_ref) = ref_target; // you can have more than one reference
ExtrudePMLLayer(mesh, axis, 0.0, ref_cond, ref_target, delta, nb_layers, false);

// instead of giving a reference for the boundary
// you can provide scal_max such that x \cdot axis = scal_max
mesh.Clear(); mesh.Read("example.mesh");

ref_cond.Zero();
Real_wp scal_max = DotProd(mesh.Vertex(0), axis); // here we compute scal_max from a vertex of the boundary
ExtrudePMLLayer(mesh, axis, scal_max, ref_cond, -1, delta, nb_layers, false);

// if the boundary is circular, axis is the center of the circle
// scal_max is the radius of the circle
// scal_max is used only if ref_target is equal to -1 (as previously)
mesh.Read("disque.mesh");

ref_cond(my_provided_ref) = ref_target; 
ExtrudePMLLayer(mesh, R2(0, 0), 0.0, ref_cond, ref_target, delta, nb_layers, true);

Location :

Mesh2D.hxx
Mesh2D.cxx

TranslateMesh

Syntax

void TranslateMesh(Mesh<Dimension2>& mesh, R2 translat)

This method translates the 2-D mesh with the vector translat.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// defining the vector of translation
R2 translat(2.5, -0.4);

// the mesh is translated
TranslateMesh(mesh, translat);

Location :

FunctionsMesh.hxx
FunctionsMesh.cxx

RotateMesh

Syntax

void RotateMesh(Mesh<Dimension2>& mesh, R2 center, Real_wp teta)

This method rotates the 2-D mesh around the point center with an angle teta.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// defining the center of rotation
R2 center(2.5, -0.4);

// the mesh is rotated by an angle pi/4
RotateMesh(mesh, center, pi_wp/4);

Location :

FunctionsMesh.hxx
FunctionsMesh.cxx

ScaleMesh

Syntax

void ScaleMesh(Mesh<Dimension2>& mesh, R2 scale)

This method scales the mesh. Each coordinate is multiplied by a given factor.

Example :

Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// defining the scaling factors for each direction
R2 scale(2.0, 3.0);

// the mesh is scaled (here x is multiplied by 2, y by 3)
ScaleMesh(mesh, scale);

Location :

FunctionsMesh.hxx
FunctionsMesh.cxx

WriteElementMesh

Syntax

void WriteElementMesh(Mesh<Dimension2>& mesh, ElementGeomReference<Dimension2>& Fb, SetPoints<Dimension2>& PointsElem, string file_name, int num_elem)

This method writes the element num_elem of the mesh in a file. The element will be subdivided depending on the geometry order. This function is useful to check a curved element.

Example :

Mesh<Dimension2> mesh;
// you can read a curved mesh
mesh.Read("test.msh");

// writes the element 4 in a new mesh (the element will be subdivided if the geometry order is greater or equal to 2).
int num_elem = 4;
const ElementGeomReference<Dimension2>& Fb = mesh.GetReferenceElement(num_elem);
SetPoints<Dimension2> PointsElem;
VectR2 s; mesh.GetVerticesElement(num_elem, s);
Fb.FjElemNodal(s, PointsElem, mesh, num_elem);
WriteElementMesh(mesh, Fb, PointsElem, "elem.mesh", 4);

Location :

FunctionsMesh.hxx
FunctionsMesh.cxx

CopyInputData

Syntax

void CopyInputData(MeshNumbering<Dimension2>& mesh_num)

This method copies the input datas of another mesh numbering. The input datas are the parameters provided in the data file.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

MeshNumberingh<Dimension2> mesh_num2(mesh);

// we copy all input data of the numbering
mesh_num2.CopyInputData(mesh_num);

Location :

NumberMesh.hxx
NumberMesh.cxx

ReallocateNeighborElements

Syntax

void ReallocateNeighborElements(int N)

This method allocates the arrays storing informations about neighbor elements.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// if you know the number of neighbor elements (elements on another processor)
mesh_num.ReallocateNeighborElements(10);

// filling informations : order, edge number, type of element
mesh_num.SetOrderNeighborElement(0, 4);
mesh_num.SetEdgeNeighborElement(0, 11);
mesh_num.SetTypeNeighborElement(0, 0);
mesh_num.SetOrderNeighborElement(1, 5);
mesh_num.SetEdgeNeighborElement(1, 27);
mesh_num.SetTypeNeighborElement(1, 0);
// ..

// once you filled informations, you call FinalizeNeighborElements()
mesh_num.FinalizeNeighborElements(mesh.GetNbBoundaryRef());

Location :

NumberMesh.hxx
NumberMesh.cxx

FinalizeNeighborElements

Syntax

void FinalizeNeighborElements(int N)

This method finalizes the arrays storing informations about neighbor elements. N is the number of referenced edges.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// if you know the number of neighbor elements (elements on another processor)
mesh_num.ReallocateNeighborElements(10);

// filling informations : order, edge number, type of element
mesh_num.SetOrderNeighborElement(0, 4);
mesh_num.SetEdgeNeighborElement(0, 11);
mesh_num.SetTypeNeighborElement(0, 0);
mesh_num.SetOrderNeighborElement(1, 5);
mesh_num.SetEdgeNeighborElement(1, 27);
mesh_num.SetTypeNeighborElement(1, 0);
// ..

// once you filled informations, you call FinalizeNeighborElements()
mesh_num.FinalizeNeighborElements(mesh.GetNbBoundaryRef());

Location :

NumberMesh.hxx
NumberMesh.cxx

GetMemorySize

Syntax

size_t GetMemorySize() const

This method returns the memory used to store the object in bytes.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

mesh_num.NumberMesh();

// size used to store the numbering ?
size_t N = mesh_num.GetMemorySize();

Location :

NumberMesh.hxx
NumberMesh.cxx

SetInputData

Syntax

void SetInputData(const string& keyword, VectString& param)

This method is used during the reading of the data file. You can look at the section devoted to the data file and see the fields related to the mesh numbering.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
// some parameters are read on a data file
ReadInputFile("test.ini", &mesh_num);

Related topics :

Data files

Location :

NumberMesh.hxx
NumberMesh.cxx

Clear

Syntax

void Clear()

This method clears the mesh numbering.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);

// we read mesh
mesh.Read("test.mesh");
// we number the mesh
mesh_num.SetOrder(2);
mesh_num.NumberMesh();

// then you can clear the mesh numbering
mesh_num.Clear();

Location :

NumberMesh.hxx
NumberMesh.cxx

Element

Syntax

ElementNumbering& Element(int i)

This method gives access to an element of the mesh numbering.

Example :


// mesh is read  
Mesh<Dimension2> mesh;
mesh.Read("test.mesh");

// mesh is numbered
MeshNumbering<Dimension2> mesh_num(mesh);

mesh_num.SetOrder(r);
mesh_num.number_map.SetOrder(mesh, r);
mesh_num.NumberMesh();

// you can retrieve dof numbers
int i = 15, j = 1;
int num_dof = mesh_num.Element(i).GetNumberDof(j);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetNbNeighborElt

Syntax

int GetNbNeighborElt() const

This method returns the number of neighbor elements.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// assuming that the mesh has been numbered and distributed
int nb_elt_neighbor = mesh_num.GetNbNeighborElt();

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetOrderNeighborElement

Syntax

int GetOrderNeighborElement(int i) const

This method returns the order of the i-th neighbor element.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// assuming that the mesh has been numbered and distributed
// you can retrieve the order of a neighbor element
int order_neighbor = mesh_num.GetOrderNeighborElement(1);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetTypeNeighborElement

Syntax

int GetTypeNeighborElement(int i) const

This method returns the type of the i-th neighbor element. This type is equal to 0 for a triangle, 1 for a quadrilateral.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// assuming that the mesh has been numbered and distributed
// you can retrieve the order of a neighbor element
int order_neighbor = mesh_num.GetOrderNeighborElement(1);
// and type of this element
int type = mesh_num.GetTypeNeighborElement(1);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetEdgeNeighborElement

Syntax

int GetEdgeNeighborElement(int i) const
int GetLocalEdgeNumberNeighborElement(int n) const

The method GetEdgeNeighborElement returns the edge number of the i-th neighbor element. The neighbor element is located on the other side of this edge. The method GetLocalEdgeNumberNeighborElement is the reciprocal of method GetEdgeNeighborElement since it returns the neighbor element number i from the edge number.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// assuming that the mesh has been numbered and distributed
// you can retrieve the order of a neighbor element
int order_neighbor = mesh_num.GetOrderNeighborElement(1);
// and type of this element
int type = mesh_num.GetTypeNeighborElement(1);
// and edge number where this neighbor element is located
int num_edge = mesh_num.GetEdgeNeighborElement(1);
// should give 1 here
int i = mesh_num.GetLocalEdgeNumberNeighborElement(num_edge);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetOrderNeighborElement

Syntax

void SetOrderNeighborElement(int i, int r)

This method sets the order for the i-th neighbor element.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// if you know the number of neighbor elements (elements on another processor)
mesh_num.ReallocateNeighborElements(10);

// filling informations : order, edge number, type of element
mesh_num.SetOrderNeighborElement(0, 4);
mesh_num.SetEdgeNeighborElement(0, 11);
mesh_num.SetTypeNeighborElement(0, 0);
mesh_num.SetOrderNeighborElement(1, 5);
mesh_num.SetEdgeNeighborElement(1, 27);
mesh_num.SetTypeNeighborElement(1, 0);
// ..

// once you filled informations, you call FinalizeNeighborElements()
mesh_num.FinalizeNeighborElements(mesh.GetNbBoundaryRef());

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetTypeNeighborElement

Syntax

void SetTypeNeighborElement(int i, int type)

This method sets the type for the i-th neighbor element.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// if you know the number of neighbor elements (elements on another processor)
mesh_num.ReallocateNeighborElements(10);

// filling informations : order, edge number, type of element
mesh_num.SetOrderNeighborElement(0, 4);
mesh_num.SetEdgeNeighborElement(0, 11);
mesh_num.SetTypeNeighborElement(0, 0);
mesh_num.SetOrderNeighborElement(1, 5);
mesh_num.SetEdgeNeighborElement(1, 27);
mesh_num.SetTypeNeighborElement(1, 0);
// ..

// once you filled informations, you call FinalizeNeighborElements()
mesh_num.FinalizeNeighborElements(mesh.GetNbBoundaryRef());

Location :

NumberMesh.hxx
NumberMeshInline.cxx

SetEdgeNeighborElement

Syntax

void SetEdgeNeighborElement(int i, int num)

This method sets the edge number for the i-th neighbor element.

Example :

Mesh<Dimension2> mesh;
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.SetOrder(4);

// if you know the number of neighbor elements (elements on another processor)
mesh_num.ReallocateNeighborElements(10);

// filling informations : order, edge number, type of element
mesh_num.SetOrderNeighborElement(0, 4);
mesh_num.SetEdgeNeighborElement(0, 11);
mesh_num.SetTypeNeighborElement(0, 0);
mesh_num.SetOrderNeighborElement(1, 5);
mesh_num.SetEdgeNeighborElement(1, 27);
mesh_num.SetTypeNeighborElement(1, 0);
// ..

// once you filled informations, you call FinalizeNeighborElements()
mesh_num.FinalizeNeighborElements(mesh.GetNbBoundaryRef());

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ReallocateElements

Syntax

void ReallocateElements(int n)

This method changes the number of elements present in the mesh numbering. On a regular use, this method should not be called, since NumberMesh will allocate the array correctly.

Location :

NumberMesh.hxx
NumberMeshInline.cxx

GetLocalVertexBoundary

Syntax

int GetLocalVertexBoundary(int type, int num_loc, int num_extremity)

This method returns the vertex number (in the reference element) of an extremity of an edge.

Parameters

type
type of element (0: triangle, 1: quadrangle)
num_loc
local edge number
num_extremity
which extremity of the edge (0 or 1)

Example :

  // for example, you want to know the first extremity of the third edge of the unit square
  int num = MeshNumbering<Dimension2>::GetLocalVertexBoundary(1, 2, 0);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ReallocatePeriodicDof

Syntax

void ReallocatePeriodicDof(int n)

This method changes the number of periodic dofs in the mesh numbering. On a regular use, this method should not be called, since NumberMesh will compute this number automatically (and allocate the arrays related to periodic degrees of freedom).

Location :

NumberMesh.hxx
NumberMesh.cxx

GetReferenceQuadrature

Syntax

const VectReal_wp& GetReferenceQuadrature(int type, int order)

This method returns the quadrature points to use for the unit interval.

Parameters

type
type of element (0: triangle, 1: quadrangle), unused here
order
order of approximation

Example :

// assuming that you constructed a mesh numbering
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// you can retrieve the quadrature points for the unit interval [0, 1]
int order = 3; // order present in the mesh
VectReal_wp points = mesh_num.GetReferenceQuadrature(0, order);

Location :

NumberMesh.hxx
NumberMeshInline.cxx

ReconstructOffsetDofs

Syntax

void ReconstructOffsetDofs()

This method reconstructs the arrays OffsetDofFaceNumber, OffsetDofEdgeNumber, OffsetDofElementNumber, OffsetDofVertexNumber. On a regular use, this method does not need to be called, since NumberMesh will call it automatically. However, if you modified dof numbers manually, it may be needed to call this method.

Example :

// assuming that you constructed a mesh numbering
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// you modify dofs manually
mesh_num.Element(1).SetNumberDof(10, 51);
// etc

// you need to reconstruct offset arrays
mesh_num.ReconstructOffsetDofs();

Location :

NumberMesh.hxx
NumberMesh.cxx

GetNbDofBoundary

Syntax

int GetNbDofBoundary(int i)

This method returns the number of degrees of freedom associated with an edge of the mesh (including degrees of freedom associated with vertices of the edge).

Example :

// assuming that you constructed a mesh numbering
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// you want to know the number of dofs on an edge
int num_edge = 42;
int nb_dof = mesh_num.GetNbDofBoundary(num_edge);

Location :

NumberMesh.hxx
NumberMesh.cxx

GetDofBoundary

Syntax

void GetDofBoundary(int i, IVect& num)

This method fills the array num with dof numbers of an edge of the mesh (including degrees of freedom associated with vertices of the edge).

Example :

// assuming that you constructed a mesh numbering
MeshNumbering<Dimension2> mesh_num(mesh);
mesh_num.NumberMesh();

// you want to know the dof numbers on an edge
int num_edge = 42; IVect num;
mesh_num.GetDofBoundary(num_edge, num);

Location :

NumberMesh.hxx
NumberMesh.cxx

SendMeshToProcessor

Syntax

void SendMeshToProcessor (Mesh<Dimension2>& glob_mesh, MeshNumbering<Dimension2>& glob_mesh_num,
IVect& NumElement, IVect& Epart, IVect& NumLoc, Mesh<Dimension2>& mesh,
int proc, MPI_Comm comm, ParamParallelMesh<Dimension2>& param, int tag)

This method sends a mesh (and its numbering) to another processor. It is intended for parallel computation where an initial global mesh is split into sub-meshes (each sub-mesh being treated on a processor).

Parameters

glob_mesh (in)
global mesh
glob_mesh_num (in)
global mesh numbering
NumElement (in)
extracted element numbers
Epart (in)
Epart(i) stores the processor number for the element i
NumLoc (in)
NumLoc(i) stores the local element number for the element i
mesh (in)
extracted mesh to send
proc (in)
processor number
comm (in)
MPI communicator
param (inout)
additional parameters
tag (optional)
tag number used in MPI communications

Example :

int nb_proc = MPI::COMM_WORLD.Get_size();
Vector<IVect> NeighboringConnectivity, SharedDofs;
IVect ProcSharedDofs;
if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    // root processors constructs a global mesh
    Mesh<Dimension2> mesh;
    mesh.Read("test.mesh");
  
    // the initial mesh is split, for example with Metis
    IVect weight_elt(mesh.GetNbElt()); weight_elt.Fill(1);
    Vector<IVect> NumElement_Subdomain;
    Vector<Mesh<Dimension2> > mesh_subdomain;
    mesh.SplitMetis(nb_proc, weight_elt, NumElement_Subdomain, mesh_subdomain, 0, false);

    // filling arrays Epart and NumLoc
    IVect Epart(mesh.GetNbElt());
    IVect NumLoc(mesh.GetNbElt());
    for (int i = 0; i < NumElement_Subdomain.GetM(); i++)
       for (int j = 0; j < NumElement_Subdomain(i).GetM(); j++)
          {
            int ne = NumElement_Subdomain(i)(j);
            Epart(ne) = i;
            NumLoc(ne) = j;
          }

    // assuming that you constructed a global mesh numbering
    MeshNumbering<Dimension2> mesh_num(mesh);
    mesh_num.NumberMesh();

    // periodic references are updated (if present)
    mesh.UpdatePeriodicReferenceSplit(Epart, NumLoc, mesh_subdomain);

    // broadcasting the global number of dofs
    int nodl_mesh = mesh_num.GetNbDof();
    MPI_Bcast(&nodl_mesh, 1, MPI_INTEGER, 0, MPI_COMM_WORLD);

    // class storing various parameters
    ParamParallelMesh<Dimension2> param;
    IVect& MinimalProc = param.MinimalProc;
    param.dg_form = ElementReference_Base::CONTINUOUS;
    param.nodl_mesh = nodl_mesh;

    // MinimalProc(i) is the minimal proc number of the dof i
    MinimalProc.Reallocate(nodl_mesh);
    MinimalProc.Fill(nb_proc+1);
    for (int i = 0; i < mesh.GetNbElt(); i++)
       for (int j = 0; j < mesh_num.Element(i).GetNbDof(); j++)
           {
              int num_dof = mesh_num.Element(i).GetNumberDof(j);
              if (num_dof >= 0)
                 MinimalProc(num_dof) = min(MinimalProc(num_dof), Epart(i)+1);
           }

     // for the root processor, the mesh is initialized
     Mesh<Dimension2> mesh_root, sub_mesh;
     MeshNumbering<Dimension2> mesh_root_num(mesh_root);
     mesh_root.CopyInputData(mesh);
     mesh_root_num.SetOrder(mesh_num.GetOrder());
     mesh_root_num.compute_dof_pml = mesh_num.compute_dof_pml;

     // loop over sub-meshes
     for (int i = 0; i < nb_proc; i++)
        {
           SendMeshToProcessor(mesh, mesh_num, NumElement_Subdomain(i), Epart, NumLoc,
                               mesh_subdomain(i), i, MPI_COMM_WORLD, param);
           if (i == 0)
             {
               mesh_root_num.number_map = mesh_num.number_map;
               RecvMeshFromProcessor(mesh_root, mesh_root_num, NeighboringConnectivity,
                                     ProcSharedDofs,  SharedDofs,
                                     param, 0, MPI_COMM_WORLD);
              }
        }

      // you can clear the global mesh
      mesh.Clear(); mesh_num.Clear();

      // for the root processor, we copy mesh_root
      mesh = mesh_root;
      mesh_num = mesh_root_num;
  }
else
  {
    // retrieving the global number of dofs
    int nodl_mesh;
    MPI_Bcast(&nodl_mesh, 1, MPI_INTEGER, 0, MPI_COMM_WORLD);

    ParamParallelMesh<Dimension2> param;
    param.nodl_mesh = nodl_mesh;
    param.dg_form = ElementReference_Base::CONTINUOUS;;
        
    // receiving the mesh associated with the processor
    mesh.ClearPeriodicCondition();
    RecvMeshFromProcessor(mesh, mesh_num, NeighboringConnectivity,
                          ProcSharedDofs,  SharedDofs, param, 0, MPI_COMM_WORLD);
    
  }

Location :

ParallelMeshFunctions.hxx
ParallelMeshFunctions.cxx

RecvMeshFromProcessor

Syntax

void RecvMeshToProcessor (Mesh<Dimension2>& mesh, MeshNumbering<Dimension2>& mesh_num,
Vector<IVect>& ConnecEdge, IVect& ProcSharedDofs, Vector<IVect>& SharedDofs,
ParamParallelMesh<Dimension2>& param, int proc, MPI_Comm comm,int tag)

This method receives a mesh (and its numbering) from another processor. It is intended for parallel computation where an initial global mesh is split into sub-meshes (each sub-mesh being treated on a processor). ConnecEdge stores various informations about neighboring elements (located on another processor). ProcSharedDofs lists all processors that share degrees of freedom with the current processor. SharedDofs(i) is the list of degrees of freedom that are shared with the processor ProcSharedDofs(i). These arrays are useful if you want to assemble a distributed matrix.

Parameters

mesh (out)
received mesh
mesh_num (out)
received mesh numbering
ConnecEdge (out)
informations about neighbor elements
ProcSharedDofs (out)
processors that share dofs with current processor
SharedDofs (out)
list of degrees of freedom shared with the processor ProcSharedDofs(i)
param (out)
class storing arrays used for MPI communications
proc (in)
processor number
comm (in)
MPI communicator

Example :

int nb_proc = MPI::COMM_WORLD.Get_size();
Vector<IVect> NeighboringConnectivity, SharedDofs;
IVect ProcSharedDofs;
if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    // root processors constructs a global mesh
    Mesh<Dimension2> mesh;
    mesh.Read("test.mesh");
  
    // the initial mesh is split, for example with Metis
    IVect weight_elt(mesh.GetNbElt()); weight_elt.Fill(1);
    Vector<IVect> NumElement_Subdomain;
    Vector<Mesh<Dimension2> > mesh_subdomain;
    mesh.SplitMetis(nb_proc, weight_elt, NumElement_Subdomain, mesh_subdomain,
0, false);

    // filling arrays Epart and NumLoc
    IVect Epart(mesh.GetNbElt());
    IVect NumLoc(mesh.GetNbElt());
    for (int i = 0; i < NumElement_Subdomain.GetM(); i++)
       for (int j = 0; j < NumElement_Subdomain(i).GetM(); j++)
          {
            int ne = NumElement_Subdomain(i)(j);
            Epart(ne) = i;
            NumLoc(ne) = j;
          }

    // assuming that you constructed a global mesh numbering
    MeshNumbering<Dimension2> mesh_num(mesh);
    mesh_num.NumberMesh();

    // periodic references are updated (if present)
    mesh.UpdatePeriodicReferenceSplit(Epart, NumLoc, mesh_subdomain);

    // broadcasting the global number of dofs
    int nodl_mesh = mesh_num.GetNbDof();
    MPI_Bcast(&nodl_mesh, 1, MPI_INTEGER, 0, MPI_COMM_WORLD);

    // class storing various parameters
    ParamParallelMesh<Dimension2> param;
    IVect& MinimalProc = param.MinimalProc;
    param.dg_form = ElementReference_Base::CONTINUOUS;
    param.nodl_mesh = nodl_mesh;

    // MinimalProc(i) is the minimal proc number of the dof i
    MinimalProc.Reallocate(nodl_mesh);
    MinimalProc.Fill(nb_proc+1);
    for (int i = 0; i < mesh.GetNbElt(); i++)
       for (int j = 0; j < mesh_num.Element(i).GetNbDof(); j++)
           {
              int num_dof = mesh_num.Element(i).GetNumberDof(j);
              if (num_dof >= 0)
                 MinimalProc(num_dof) = min(MinimalProc(num_dof), Epart(i)+1);
           }

     // for the root processor, the mesh is initialized
     Mesh<Dimension2> mesh_root, sub_mesh;
     MeshNumbering<Dimension2> mesh_root_num(mesh_root);
     mesh_root.CopyInputData(mesh);
     mesh_root_num.SetOrder(mesh_num.GetOrder());
     mesh_root_num.compute_dof_pml = mesh_num.compute_dof_pml;

     // loop over sub-meshes
     for (int i = 0; i < nb_proc; i++)
        {
           SendMeshToProcessor(mesh, mesh_num, NumElement_Subdomain(i), Epart, NumLoc,
                               mesh_subdomain(i), i, MPI_COMM_WORLD, param);
           if (i == 0)
             {
               mesh_root_num.number_map = mesh_num.number_map;
               RecvMeshFromProcessor(mesh_root, mesh_root_num, NeighboringConnectivity,
                                     ProcSharedDofs,  SharedDofs,
                                     param, 0, MPI_COMM_WORLD);
              }
         }

      // you can clear the global mesh
      mesh.Clear(); mesh_num.Clear();

      // for the root processor, we copy mesh_root
      mesh = mesh_root;
      mesh_num = mesh_root_num;
  }
else
  {
    // retrieving the global number of dofs
    int nodl_mesh;
    MPI_Bcast(&nodl_mesh, 1, MPI_INTEGER, 0, MPI_COMM_WORLD);

    ParamParallelMesh<Dimension2> param;
    param.nodl_mesh = nodl_mesh;
    param.dg_form = ElementReference_Base::CONTINUOUS;;
        
    // receiving the mesh associated with the processor
    mesh.ClearPeriodicCondition();
    RecvMeshFromProcessor(mesh, mesh_num, NeighboringConnectivity,
                          ProcSharedDofs,  SharedDofs, param, 0, MPI_COMM_WORLD);

  }

Location :

ParallelMeshFunctions.hxx
ParallelMeshFunctions.cxx