|
NGSolve
4.9
|
Access to mesh topology and geometry. More...
#include <meshaccess.hpp>
Public Member Functions | |
| MeshAccess () | |
| connects to Netgen - mesh | |
| virtual | ~MeshAccess () |
| not much to do | |
| int | GetDimension () const |
| the spatial dimension of the mesh | |
| int | GetNP () const |
| number of points. needed for isoparametric nodal elements | |
| int | GetNV () const |
| number of vertices | |
| int | GetNE () const |
| number of elements in the domain | |
| int | GetNSE () const |
| number of boundary elements | |
| int | GetNEdges () const |
| number of edges in the whole mesh | |
| int | GetNFaces () const |
| number of faces in the whole mesh | |
| int | GetNDomains () const |
| maximal sub-domain (material) index. range is [0, ndomains) | |
| int | GetNBoundaries () const |
| maximal boundary condition index. range is [0, nboundaries) | |
| template<int D> | |
| void | GetPoint (int pi, Vec< D > &p) const |
| returns point coordinate | |
| template<int D> | |
| Vec< D > | GetPoint (int pi) const |
| returns point coordinate | |
| ELEMENT_TYPE | GetElType (int elnr) const |
| the geometry type of the element | |
| ELEMENT_TYPE | GetSElType (int elnr) const |
| the geometry type of the boundary element | |
| int | GetElIndex (int elnr) const |
| the sub-domain index of the element | |
| int | GetSElIndex (const int elnr) const |
| the boundary-condition index of the boundary element | |
| void | SetElIndex (int elnr, int index) const |
| change sub-domain index (???) | |
| string | GetElMaterial (int elnr) const |
| the material of the element | |
| string | GetDomainMaterial (int domnr) const |
| the material of the sub-domain | |
| string | GetSElBCName (int selnr) const |
| the boundary condition name of surface element | |
| string | GetBCNumBCName (const int bcnr) const |
| the boundary condition name of boundary condition number | |
| int | GetSElSurfaceNumber (const int elnr) const |
| not sure who needs that | |
| int | GetSElFDNumber (const int elnr) const |
| not sure who needs that | |
| void | GetSElNeighbouringDomains (const int elnr, int &in, int &out) const |
| the sub-domain indices next to boundary element. | |
| void | UpdateBuffers () |
| update buffered global quantities. | |
| int | GetNElements (int dim) const |
| number of elements of dimension dim | |
| int | GetNNodes (NODE_TYPE nt) const |
| number of nodes of type nt | |
| void | GetTopologicElement (int elnr, TopologicElement &topel) const |
| the topology of a domain - element experimental, not recommended for use yet | |
| Ng_Element | GetElement (int elnr) const |
| returns topology of a Netgen - element. | |
| Ng_Element | GetSElement (int elnr) const |
| returns topology of a Netgen - element. | |
| template<int DIM> | |
| Ng_Element | GetElement (int elnr) const |
| returns element of compile-time fixed dimension | |
| template<int DIM> | |
| Ng_Node< DIM > | GetNode (int nr) const |
| returns node topology. | |
| int | GetElNV (int elnr) const |
| void | GetElPNums (int elnr, Array< int > &pnums) const |
| returns the points of an element. | |
| void | GetSElPNums (int selnr, Array< int > &pnums) const |
| returns the points of a boundary element. | |
| void | GetElVertices (int elnr, Array< int > &vnums) const |
| returns the vertices of an element | |
| void | GetSElVertices (int selnr, Array< int > &vnums) const |
| returns the vertices of a boundary element | |
| void | GetElEdges (int elnr, Array< int > &ednums) const |
| returns the edges of an element | |
| void | GetElEdges (int elnr, Array< int > &ednums, Array< int > &orient) const |
| void | GetSElEdges (int selnr, Array< int > &ednums) const |
| returns the edges of a boundary element | |
| void | GetSElEdges (int selnr, Array< int > &ednums, Array< int > &orient) const |
| void | GetElFaces (int elnr, Array< int > &fnums) const |
| returns the faces of an element | |
| void | GetElFaces (int elnr, Array< int > &fnums, Array< int > &orient) const |
| int | GetSElFace (int selnr) const |
| returns face number of surface element | |
| void | GetSElFace (int selnr, int &fnum, int &orient) const |
| void | GetFacePNums (int fnr, Array< int > &pnums) const |
| returns vertex numbers of face | |
| void | GetEdgePNums (int enr, int &pn1, int &pn2) const |
| returns vertex numbers of edge | |
| void | GetEdgePNums (int enr, Array< int > &pnums) const |
| returns vertex numbers of edge | |
| void | GetEdgeElements (int enr, Array< int > &elnums) const |
| returns all elements connected to an edge | |
| void | GetFaceEdges (int fnr, Array< int > &edges) const |
| returns all edges of a face | |
| void | GetFaceElements (int fnr, Array< int > &elnums) const |
| returns elements connected to a face | |
| void | GetSegmentPNums (int snr, Array< int > &pnums) const |
| point numbers of a 1D element | |
| int | GetSegmentIndex (int snr) const |
| index of 1D element | |
| int | GetNFacets () const |
| number of facets of an element. | |
| void | GetElFacets (int elnr, Array< int > &fnums) const |
| facets of an element | |
| void | GetSElFacets (int selnr, Array< int > &fnums) const |
| facet of a surface element | |
| void | GetFacetPNums (int fnr, Array< int > &pnums) const |
| vertices of a facet | |
| ELEMENT_TYPE | GetFacetType (int fnr) const |
| geometry type of facet | |
| void | GetFacetElements (int fnr, Array< int > &elnums) const |
| elements connected to facet | |
| int | GetElOrder (int enr) const |
| element order stored in Netgen | |
| INT< 3 > | GetElOrders (int enr) const |
| anisotropic order stored in Netgen | |
| void | SetElOrder (int enr, int order) const |
| set element order in Netgen | |
| void | SetElOrders (int enr, int ox, int oy, int oz) const |
| set anisotropic element order in Netgen | |
| int | GetSElOrder (int enr) const |
| order of suface element | |
| INT< 2 > | GetSElOrders (int enr) const |
| anisotropic order of suface element | |
| void | SetSElOrder (int enr, int order) const |
| set surface element order | |
| void | SetSElOrders (int enr, int ox, int oy) const |
| set anisotropic surface element order | |
| double | ElementVolume (int elnr) const |
| the volume of an element (mid-point rule) | |
| double | SurfaceElementVolume (int selnr) const |
| the area of a boundary element (mid-point rule) | |
| int | GetNLevels () const |
| number of multigrid levels | |
| void | GetParentNodes (int pi, int *parents) const |
| the two parent vertices of a vertex. -1 for coarse-grid vertices | |
| int | GetParentElement (int ei) const |
| number of parent element on next coarser mesh | |
| int | GetParentSElement (int ei) const |
| number of parent boundary element on next coarser mesh | |
| int | GetClusterRepVertex (int pi) const |
| representant of vertex for anisotropic meshes | |
| int | GetClusterRepEdge (int pi) const |
| representant of edge for anisotropic meshes | |
| int | GetClusterRepFace (int pi) const |
| representant of face for anisotropic meshes | |
| int | GetClusterRepElement (int pi) const |
| representant of element for anisotropic meshes | |
| ngfem::ElementTransformation & | GetTrafo (int elnr, bool boundary, LocalHeap &lh) const |
| returns the transformation from the reference element to physical element. | |
| void | SetPointSearchStartElement (const int el) const |
| int | FindElementOfPoint (FlatVector< double > point, IntegrationPoint &ip, bool build_searchtree, const Array< int > *const indices=NULL) const |
| int | FindElementOfPoint (FlatVector< double > point, ngfem::IntegrationPoint &ip, bool build_searchtree, const int index) const |
| int | FindSurfaceElementOfPoint (FlatVector< double > point, ngfem::IntegrationPoint &ip, bool build_searchtree, const Array< int > *const indices=NULL) const |
| int | FindSurfaceElementOfPoint (FlatVector< double > point, ngfem::IntegrationPoint &ip, bool build_searchtree, const int index) const |
| bool | IsElementCurved (int elnr) const |
| is element straiht or curved ? | |
| void | GetPeriodicVertices (Array< ngstd::INT< 2 > > &pairs) const |
| int | GetNPairsPeriodicVertices () const |
| void | GetPeriodicVertices (int idnr, Array< ngstd::INT< 2 > > &pairs) const |
| int | GetNPairsPeriodicVertices (int idnr) const |
| void | GetPeriodicEdges (Array< ngstd::INT< 2 > > &pairs) const |
| int | GetNPairsPeriodicEdges () const |
| void | GetPeriodicEdges (int idnr, Array< ngstd::INT< 2 > > &pairs) const |
| int | GetNPairsPeriodicEdges (int idnr) const |
| virtual void | PushStatus (const char *str) const |
| virtual void | PopStatus () const |
| virtual void | SetThreadPercentage (double percent) const |
| virtual void | GetStatus (string &str, double &percent) const |
| virtual void | SetTerminate (void) const |
| virtual void | UnSetTerminate (void) const |
| virtual bool | ShouldTerminate (void) const |
| void | GetVertexElements (int vnr, Array< int > &elems) const |
| void | GetVertexSurfaceElements (int vnr, Array< int > &elems) const |
| void | SetHigherIntegrationOrder (int elnr) |
| void | UnSetHigherIntegrationOrder (int elnr) |
| void | LoadMeshFromString (const string &str) |
| void | InitPointCurve (double red=1, double green=0, double blue=0) const |
| void | AddPointCurvePoint (const Vec< 3 > &point) const |
| MPI_Comm | GetCommunicator () const |
| void | GetDistantProcs (Node node, Array< int > &procs) const |
| Returns the list of other MPI - processes where node is present. | |
| int | GetGlobalNodeNum (Node node) const |
| Returns the global number of the node. | |
| template<typename T > | |
| void | AllReduceNodalData (NODE_TYPE nt, Array< T > &data, MPI_Op op) const |
| Reduces MPI - distributed data associated with mesh-nodes. | |
Access to mesh topology and geometry.
MeshAccess provides information such as element-to-vertex table, elemenet-to-edge table etc.
It provides also element-geometry (straight or curved elements) via GetTrafo.
Internally, MeshAccess calls functions from Netgen.
| void ngcomp::MeshAccess::GetDistantProcs | ( | Node | node, |
| Array< int > & | procs | ||
| ) | const |
Returns the list of other MPI - processes where node is present.
The ordering coincides for all processes.
| Ng_Element ngcomp::MeshAccess::GetElement | ( | int | elnr | ) | const [inline] |
returns topology of a Netgen - element.
This is the new (2008), unified concept. The function returns a direct access to the Netgen mesh structure instead of copying point numbers etc. The nasty 1-0 convertion is done on the fly.
| void ngcomp::MeshAccess::GetElPNums | ( | int | elnr, |
| Array< int > & | pnums | ||
| ) | const |
returns the points of an element.
vertices and possibly edge-midpoints
| int ngcomp::MeshAccess::GetGlobalNodeNum | ( | Node | node | ) | const |
Returns the global number of the node.
Currently, this function works only for vertex-nodes.
| int ngcomp::MeshAccess::GetNFacets | ( | ) | const [inline] |
number of facets of an element.
facets are edges (2D) or faces (3D)
| Ng_Node<DIM> ngcomp::MeshAccess::GetNode | ( | int | nr | ) | const [inline] |
returns node topology.
A facet or edge-node knows its vertices etc. The method is not yet fully functional.
| Ng_Element ngcomp::MeshAccess::GetSElement | ( | int | elnr | ) | const [inline] |
returns topology of a Netgen - element.
This is the new (2008), unified concept. The function returns a direct access to the Netgen mesh structure instead of copying point numbers etc. The nasty 1-0 convertion is done on the fly.
| void ngcomp::MeshAccess::GetSElNeighbouringDomains | ( | const int | elnr, |
| int & | in, | ||
| int & | out | ||
| ) | const [inline] |
the sub-domain indices next to boundary element.
returns -1 for void
| void ngcomp::MeshAccess::GetSElPNums | ( | int | selnr, |
| Array< int > & | pnums | ||
| ) | const |
returns the points of a boundary element.
vertex and possibly edge-midpoints
| ngfem::ElementTransformation& ngcomp::MeshAccess::GetTrafo | ( | int | elnr, |
| bool | boundary, | ||
| LocalHeap & | lh | ||
| ) | const |
returns the transformation from the reference element to physical element.
Given a point in the refrence element, the ElementTransformation can compute the mapped point as well as the Jacobian
| void ngcomp::MeshAccess::UpdateBuffers | ( | ) |
update buffered global quantities.
Must be called after every change of the mesh topology
1.7.6.1