Actual source code: meshQuery.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: meshQuery.c,v 1.19 2000/10/17 13:48:57 knepley Exp $";
3: #endif
5: /*
6: Defines the interface to the mesh query functions
7: */
9: #include src/mesh/meshimpl.h
11: /*----------------------------------------------- Mesh Query Functions -----------------------------------------------*/
12: /*@
13: MeshSetDimension - This function sets the dimension of the Mesh. This may only be done before a mesh is generated.
15: Not collective
17: Input Parameters:
18: + mesh - The mesh
19: - dim - The mesh dimension
21: Level: intermediate
23: .keywords: Mesh, dimension
24: .seealso: MatGetDimension(), PartitionGetDimension()
25: @*/
26: int MeshSetDimension(Mesh mesh, int dim)
27: {
30: if ((dim < 0) || (dim > 3)) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Dimension %d is invalid", dim);
31: if (mesh->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Dimension cannot be set after mesh has been generated");
32: mesh->dim = dim;
33: return(0);
34: }
36: /*@
37: MeshGetDimension - This function returns the dimension of the Mesh
39: Not collective
41: Input Parameter:
42: . mesh - The mesh
44: Output Parameter:
45: . dim - The mesh dimension
47: Level: intermediate
49: .keywords: Mesh, dimension
50: .seealso: MatSetDimension(), PartitionGetDimension()
51: @*/
52: int MeshGetDimension(Mesh mesh, int *dim)
53: {
57: *dim = mesh->dim;
58: return(0);
59: }
61: /*@
62: MeshGetInfo - This function returns some information about the mesh.
64: Collective on Mesh
66: Input Parameter:
67: . mesh - The mesh
69: Output Parameters:
70: + vertices - The number of vertices
71: . nodes - The number of nodes
72: . edges - The number of edges
73: - elements - The number of elements
75: Level: intermediate
77: .keywords: mesh
78: .seealso: MatGetCorners()
79: @*/
80: int MeshGetInfo(Mesh mesh, int *vertices, int *nodes, int *edges, int *elements) {
83: if (vertices != PETSC_NULL) *vertices = mesh->numVertices;
84: if (nodes != PETSC_NULL) *nodes = mesh->numNodes;
85: if (edges != PETSC_NULL) *edges = mesh->numEdges;
86: if (elements != PETSC_NULL) *elements = mesh->numFaces;
87: return(0);
88: }
90: /*@
91: MeshSetNumCorners - This function sets the number of nodes on each element. This may only be done before a mesh is generated.
93: Not collective
95: Input Parameters:
96: + mesh - The mesh
97: - numCorners - The number of nodes on each element
99: Level: intermediate
101: .keywords: Mesh, dimension
102: .seealso: MatGetNumCorners()
103: @*/
104: int MeshSetNumCorners(Mesh mesh, int numCorners)
105: {
108: if (numCorners <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Number of corners %d is invalid", numCorners);
109: if (mesh->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Corners cannot be set after mesh has been generated");
110: mesh->numCorners = numCorners;
111: return(0);
112: }
114: /*@
115: MeshGetNumCorners - This function returns the number of nodes on each element.
117: Not collective
119: Input Parameter:
120: . mesh - The mesh
122: Output Parameter:
123: . numCorners - The number of nodes on an element
125: Level: intermediate
127: .keywords: mesh, corner
128: .seealso: MatGetInfo()
129: @*/
130: int MeshGetNumCorners(Mesh mesh, int *numCorners) {
134: *numCorners = mesh->numCorners;
135: return(0);
136: }
138: /*@
139: MeshSetBoundingBox - This function sets the bounding box for the mesh. This can only be done before
140: the mesh is generated.
142: Not collective
144: Input Parameters:
145: + mesh - The mesh
146: . startX, startY, startZ - The lower-left corner of the box
147: - endX, endY, endZ - The upper-right corner of the box
149: Level: intermediate
151: .keywords: mesh, bounding box
152: .seealso: MeshGetBoundingBox(), MeshSetLocalBoundingBox(), MeshUpdateBoundingBox()
153: @*/
154: int MeshSetBoundingBox(Mesh mesh, PetscReal startX, PetscReal startY, PetscReal startZ, PetscReal endX, PetscReal endY, PetscReal endZ)
155: {
158: if ((startX > endX) || (startY > endY) || (startZ > endZ)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "The bounding box is inverted");
159: if (mesh->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Bounding box cannot be set after mesh has been generated");
160: mesh->startX = startX;
161: mesh->startY = startY;
162: mesh->startZ = startZ;
163: mesh->endX = endX;
164: mesh->endY = endY;
165: mesh->endZ = endZ;
166: mesh->sizeX = endX - startX;
167: mesh->sizeY = endY - startY;
168: mesh->sizeZ = endZ - startZ;
169: return(0);
170: }
172: /*@
173: MeshGetBoundingBox - This function returns the bounding box for the mesh
175: Not collective
177: Input Parameter:
178: . mesh - The mesh
180: Output Parameters:
181: + startX, startY, startZ - The lower-left corner of the box
182: - endX, endY, endZ - The upper-right corner of the box
184: Level: intermediate
186: .keywords: mesh, bounding box
187: .seealso: MeshSetBoundingBox(), MeshGetLocalBoundingBox(), MeshUpdateBoundingBox()
188: @*/
189: int MeshGetBoundingBox(Mesh mesh, PetscReal *startX, PetscReal *startY, PetscReal *startZ, PetscReal *endX, PetscReal *endY, PetscReal *endZ)
190: {
193: if (startX != PETSC_NULL) {
195: *startX = mesh->startX;
196: }
197: if (startY != PETSC_NULL) {
199: *startY = mesh->startY;
200: }
201: if (startZ != PETSC_NULL) {
203: *startZ = mesh->startZ;
204: }
205: if (endX != PETSC_NULL) {
207: *endX = mesh->endX;
208: }
209: if (endY != PETSC_NULL) {
211: *endY = mesh->endY;
212: }
213: if (endZ != PETSC_NULL) {
215: *endZ = mesh->endZ;
216: }
217: return(0);
218: }
220: /*@
221: MeshSetLocalBoundingBox - This function sets the local bounding box for the mesh. This can only be done before
222: the mesh is generated. The local box is the smallest rectangle enclosing the portion of the mesh allocated to
223: this processor.
225: Not collective
227: Input Parameters:
228: + mesh - The mesh
229: . startX, startY, startZ - The lower-left corner of the local box
230: - endX, endY, endZ - The upper-right corner of the local box
232: Level: intermediate
234: .keywords: mesh, bounding box, local
235: .seealso: MeshGetLocalBoundingBox(), MeshSetBoundingBox(), MeshUpdateBoundingBox()
236: @*/
237: int MeshSetLocalBoundingBox(Mesh mesh, PetscReal startX, PetscReal startY, PetscReal startZ, PetscReal endX, PetscReal endY, PetscReal endZ)
238: {
241: if ((startX > endX) || (startY > endY) || (startZ > endZ)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "The bounding box is inverted");
242: if (mesh->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Bounding box cannot be set after mesh has been generated");
243: mesh->locStartX = startX;
244: mesh->locStartY = startY;
245: mesh->locStartZ = startZ;
246: mesh->locEndX = endX;
247: mesh->locEndY = endY;
248: mesh->locEndZ = endZ;
249: mesh->locSizeX = mesh->locEndX - mesh->locStartX;
250: mesh->locSizeY = mesh->locEndY - mesh->locStartY;
251: mesh->locSizeZ = mesh->locEndZ - mesh->locStartZ;
252: return(0);
253: }
255: /*@
256: MeshGetLocalBoundingBox - This function returns the local bounding box for the mesh. The local box is the smallest
257: rectangle enclosing the portion of the mesh allocated to this processor.
259: Not collective
261: Input Parameter:
262: . mesh - The mesh
264: Output Parameters:
265: + startX, startY, startZ - The lower-left corner of the local box
266: - endX, endY, endZ - The upper-right corner of the local box
268: Level: intermediate
270: .keywords: mesh, bounding box, local
271: .seealso: MeshSetLocalBoundingBox(), MeshGetBoundingBox(), MeshUpdateBoundingBox()
272: @*/
273: int MeshGetLocalBoundingBox(Mesh mesh, PetscReal *startX, PetscReal *startY, PetscReal *startZ, PetscReal *endX, PetscReal *endY, PetscReal *endZ)
274: {
277: if (startX != PETSC_NULL) {
279: *startX = mesh->locStartX;
280: }
281: if (startY != PETSC_NULL) {
283: *startY = mesh->locStartY;
284: }
285: if (startZ != PETSC_NULL) {
287: *startZ = mesh->locStartZ;
288: }
289: if (endX != PETSC_NULL) {
291: *endX = mesh->locEndX;
292: }
293: if (endY != PETSC_NULL) {
295: *endY = mesh->locEndY;
296: }
297: if (endZ != PETSC_NULL) {
299: *endZ = mesh->locEndZ;
300: }
301: return(0);
302: }
304: /*@
305: MeshUpdateBoundingBox - This function updates the bounding box for the mesh based upon the current geometry.
307: Not collective
309: Input Parameter:
310: . mesh - The mesh
312: Level: intermediate
314: .keywords: mesh, bounding box, update
315: .seealso: MeshSetBoundingBox(), MeshGetBoundingBox()
316: @*/
317: int MeshUpdateBoundingBox(Mesh mesh)
318: {
319: MPI_Comm comm;
320: PetscTruth isPeriodic, isSetup;
321: PetscReal locStartX, locStartY, locStartZ;
322: PetscReal locEndX, locEndY, locEndZ;
323: PetscReal startX, startY, startZ;
324: PetscReal endX, endY, endZ;
325: int ierr;
329: PetscObjectGetComm((PetscObject) mesh, &comm);
330: MeshIsPeriodic(mesh, &isPeriodic);
331: if (isPeriodic == PETSC_TRUE) return(0);
332: /* Calculate the local bounding box */
333: (*mesh->ops->updateboundingbox)(mesh);
334: /* Calculate global bounding box */
335: MeshGetLocalBoundingBox(mesh, &locStartX, &locStartY, &locStartZ, &locEndX, &locEndY, &locEndZ);
336: MPI_Allreduce(&locStartX, &startX, 1, MPIU_REAL, MPI_MIN, comm);
337: MPI_Allreduce(&locEndX, &endX, 1, MPIU_REAL, MPI_MAX, comm);
338: MPI_Allreduce(&locStartY, &startY, 1, MPIU_REAL, MPI_MIN, comm);
339: MPI_Allreduce(&locEndY, &endY, 1, MPIU_REAL, MPI_MAX, comm);
340: MPI_Allreduce(&locStartZ, &startZ, 1, MPIU_REAL, MPI_MIN, comm);
341: MPI_Allreduce(&locEndZ, &endZ, 1, MPIU_REAL, MPI_MAX, comm);
342: /* Pretend that the mesh is not setup */
343: isSetup = mesh->setupcalled;
344: mesh->setupcalled = PETSC_FALSE;
345: MeshSetBoundingBox(mesh, startX, startY, startZ, endX, endY, endZ);
346: mesh->setupcalled = isSetup;
347: return(0);
348: }
350: /*@
351: MeshGetPartition - Returns the Partition object for this Mesh.
353: Not collective
355: Input Parameter:
356: . mesh - The mesh
358: Output Parameter:
359: . part - The partition
361: Level: intermediate
363: .keywords: Mesh, Partition, get
364: .seealso: PartitionGetMesh()
365: @*/
366: int MeshGetPartition(Mesh mesh, Partition *part)
367: {
371: *part = mesh->part;
372: return(0);
373: }
374: /*@
375: MeshSetMovement - This function sets the dimension of the Mesh. This may only be done before a mesh is generated.
377: Not collective
379: Input Parameters:
380: + mesh - The mesh
381: - dim - The mesh dimension
383: Level: intermediate
385: .keywords: Mesh, dimension
386: .seealso: MeshGetMovement(), MeshGetMover()
387: @*/
388: int MeshSetMovement(Mesh mesh, PetscTruth isMoving)
389: {
392: mesh->isMoving = isMoving;
393: return(0);
394: }
396: /*@
397: MeshGetMovement - This function determines whether the Mesh is moving.
399: Not collective
401: Input Parameter:
402: . mesh - The mesh
404: Output Parameter:
405: . isMoving - The flag for mesh movement
407: Level: intermediate
409: .keywords: Mesh, movement
410: .seealso: MeshSetMovement(), MeshGetMover()
411: @*/
412: int MeshGetMovement(Mesh mesh, PetscTruth *isMoving)
413: {
417: *isMoving = mesh->isMoving;
418: return(0);
419: }
421: /*------------------------------------------- Mesh Boundary Query Functions ------------------------------------------*/
422: /*@
423: MeshGetNumBoundaries - Returns the number of boundaries in the mesh.
425: Not collective
427: Input Parameter:
428: . mesh - The mesh
430: Output Parameter:
431: . numBd - The number boundaries
433: Level: intermediate
435: .keywords: mesh, boundary
436: .seealso: MeshGetBoundaryStart(), MeshGetBoundarySize()
437: @*/
438: int MeshGetNumBoundaries(Mesh mesh, int *numBd) {
442: *numBd = mesh->numBd;
443: return(0);
444: }
446: /*@
447: MeshGetBoundarySize - Returns the number of nodes on the specified boundary.
449: Not collective
451: Input Parameters:
452: + mesh - The mesh
453: - boundary - The boundary marker
455: Output Parameter:
456: . size - The number of nodes on the boundary
458: Level: intermediate
460: .keywords: mesh, boundary
461: .seealso: MeshGetBoundaryStart()
462: @*/
463: int MeshGetBoundarySize(Mesh mesh, int boundary, int *size)
464: {
470: (*mesh->ops->getboundarysize)(mesh, boundary, size);
471: return(0);
472: }
474: /*@
475: MeshGetBoundaryIndex - Returns the index of the specified boundary.
477: Not collective
479: Input Parameters:
480: + mesh - The mesh
481: - boundary - The boundary marker
483: Output Parameter:
484: . index - The index of the boundary
486: Level: intermediate
488: .keywords: mesh, boundary
489: .seealso: MeshGetBoundarySize()
490: @*/
491: int MeshGetBoundaryIndex(Mesh mesh, int boundary, int *index)
492: {
498: (*mesh->ops->getboundaryindex)(mesh, boundary, index);
499: return(0);
500: }
502: /*@
503: MeshGetBoundaryStart - Retrieves the canonical node number of the first node
504: on the specified boundary.
506: Not collective
508: Input Parameters:
509: + mesh - The mesh
510: . boundary - The boundary marker
511: - ghost - Flag for including ghost nodes
513: Output Parameter:
514: . node - The canonical node number
516: Note:
517: $ This is typically used as an iteration construct with MeshGetBoundaryNext(),
518: $ for example:
519: $
520: $ MeshGetBoundaryStart(mesh, OUTER_BD, &node, PETSC_FALSE);
521: $ while (node >= 0) {
522: $ PetscPrintf(PETSC_COMM_SELF, "I have boundary node %dn", node);
523: $ MeshGetBoundaryNext(mesh, OUTER_BD, &node, PETSC_FALSE);
524: $ }
526: Level: intermediate
528: .keywords: mesh, boundary, iterator
529: .seealso: MeshGetBoundaryNext()
530: @*/
531: int MeshGetBoundaryStart(Mesh mesh, int boundary, PetscTruth ghost, int *node)
532: {
538: (*mesh->ops->getboundarystart)(mesh, boundary, ghost, node);
539: return(0);
540: }
542: /*@
543: MeshGetBoundaryNext - Retrieves the canonical node number of the next node
544: on the specified boundary, with -1 indicating boundary termination.
546: Not collective
548: Input Parameters:
549: + mesh - The mesh
550: . boundary - The boundary marker
551: - ghost - Flag for including ghost nodes
553: Output Parameter:
554: . node - The canonical node number
556: Note:
557: $ This is typically used as an iteration construct with MeshGetBoundaryStart(),
558: $ for example:
559: $
560: $ MeshGetBoundaryStart(mesh, OUTER_BD, &node, PETSC_FALSE);
561: $ while (node >= 0) {
562: $ PetscPrintf(PETSC_COMM_SELF, "I have boundary node %dn", node);
563: $ MeshGetBoundaryNext(mesh, OUTER_BD, &node, PETSC_FALSE);
564: $ }
566: Also, this only returns nodes which lie in the given domain, thus the above
567: loop would run in parallel on all processors, returning different nodes on each.
569: Level: intermediate
571: .keywords: mesh, boundary, iterator
572: .seealso: MeshGetBoundaryStart()
573: @*/
574: int MeshGetBoundaryNext(Mesh mesh, int boundary, PetscTruth ghost, int *node)
575: {
581: (*mesh->ops->getboundarynext)(mesh, boundary, ghost, node);
582: return(0);
583: }
585: /*@
586: MeshGetActiveBoundary - Retrieves the boundary marker for the boundary
587: last iterated over, or -1 if no iteration has taken place.
589: Not collective
591: Input Parameter:
592: . mesh - The mesh
594: Output Parameter:
595: . boundary - The boundary marker
597: Level: advanced
599: .keywords: grid, boundary, iterator
600: .seealso: GridGetBoundaryNext(), MeshGetBoundaryStart()
601: @*/
602: int MeshGetActiveBoundary(Mesh mesh, int *boundary)
603: {
609: (*mesh->ops->getactiveboundary)(mesh, boundary);
610: return(0);
611: }
613: /*--------------------------------------------- Mesh Node Query Functions --------------------------------------------*/
614: /*@
615: MeshGetNodeBoundary - Returns the boundary on which the node lies, or 0 for interior nodes.
617: Not collective
619: Input Parameters:
620: + mesh - The mesh
621: - node - The node
623: Output Parameter:
624: . bd - The boundary
626: Level: intermediate
628: .keywords: mesh, boundary, node
629: .seealso: MeshGetBoundarySize()
630: @*/
631: int MeshGetNodeBoundary(Mesh mesh, int node, int *bd)
632: {
638: (*mesh->ops->getnodeboundary)(mesh, node, bd);
639: return(0);
640: }
642: /*@
643: MeshNodeIsVertex - This function determines whether or not a node is a vertex.
645: Not collective
647: Input Parameters:
648: + mesh - The mesh
649: - node - The node
651: Output Parameter:
652: . isVertex - The flag for a vertex
654: Note:
655: Vertices are nodes which are connected to edges.
657: Level: intermediate
659: .keywords: mesh, node, vertex
660: .seealso: MeshGetNodeCoords()
661: @*/
662: int MeshNodeIsVertex(Mesh mesh, int node, PetscTruth *isVertex)
663: {
669: (*mesh->ops->nodeisvertex)(mesh, node, isVertex);
670: return(0);
671: }
673: /*@
674: MeshGetNodeCoords - Retrieves the coordinates of a selected node
676: Input Parameters:
677: + mesh - The mesh
678: - node - The local node number
680: Output Parameters:
681: . x,y,z - The coordinates
683: Level: intermediate
685: .keywords: mesh, hole, coordinates
686: .seealso: MeshGetBoundaryStart()
687: @*/
688: int MeshGetNodeCoords(Mesh mesh, int node, double *x, double *y, double *z)
689: {
697: (*mesh->ops->getnodecoords)(mesh, node, x, y, z);
698: return(0);
699: }
701: /*@
702: MeshSetNodeCoords - Sets the coordinates of a selected node
704: Collective on Mesh
706: Input Parameters:
707: + mesh - The mesh
708: . node - The local node number
709: - x,y,z - The coordinates
711: Level: intermediate
713: .keywords: mesh, node, coordinates
714: .seealso: MeshGetBoundaryStart()
715: @*/
716: int MeshSetNodeCoords(Mesh mesh, int node, double x, double y, double z)
717: {
722: (*mesh->ops->setnodecoords)(mesh, node, x, y, z);
723: return(0);
724: }
726: /*@
727: MeshGetNodeCoordsSaved - Retrieves the coordinates of a selected node
728: from those saved with MeshSaveMesh().
730: Input Parameters:
731: + mesh - The mesh
732: - node - The canonical node number
734: Output Parameters:
735: . x,y,z - The coordinates
737: Level: intermediate
739: .keywords: mesh, node, coordinates
740: .seealso: MeshGetNodeCoords(), MeshSaveMesh()
741: @*/
742: int MeshGetNodeCoordsSaved(Mesh mesh, int node, double *x, double *y, double *z)
743: {
751: (*mesh->ops->getnodecoordssaved)(mesh, node, x, y, z);
752: return(0);
753: }
755: /*@
756: MeshGetNearestNode - This function returns the node nearest to
757: the given point, or -1 if the closest node is not contained in
758: the local mesh.
760: Not collective
762: Input Parameters:
763: + mesh - The mesh
764: . x,y,z - The node coordinates
765: - outside - A flag to allow points outside the domain
767: Output Parameter:
768: . node - The nearest node
770: Note:
771: The outside flag allows points outside the domain to be tested. If this flag
772: is PETSC_FALSE, and (x,y,z) does not lie in the global domain, then an error
773: will result.
775: Level: beginner
777: .keywords: mesh, node, point location
778: .seealso MeshLocatePoint()
779: @*/
780: int MeshGetNearestNode(Mesh mesh, double x, double y, double z, PetscTruth outside, int *node)
781: {
787: (*mesh->ops->nearestnode)(mesh, x, y, z, outside, node);
788: return(0);
789: }
791: /*@
792: MeshGetNearestBdNode - This function returns the boundary node nearest to
793: the given point, or -1 if the closest node is not contained in the local mesh.
795: Not collective
797: Input Parameters:
798: + mesh - The mesh
799: - x,y,z - The node coordinates
801: Output Parameter:
802: . node - The nearest boundary node
804: Level: beginner
806: .keywords: mesh, node, point location
807: .seealso MeshGetNearestNode(), MeshLocatePoint()
808: @*/
809: int MeshGetNearestBdNode(Mesh mesh, double x, double y, double z, int *node)
810: {
811: Partition p;
812: double minDist = PETSC_MAX;
813: double dist, nx, ny, nz;
814: int minNode, globalMinNode, allMinNode;
815: int numCorners, startNode, endNode;
816: int elem, corner, nearNode, bd;
817: int ierr;
822: MeshGetNumCorners(mesh, &numCorners);
823: MeshGetPartition(mesh, &p);
824: PartitionGetStartNode(p, &startNode);
825: PartitionGetEndNode(p, &endNode);
826: MeshLocatePoint(mesh, x, y, z, &elem);
827: if (elem >= 0) {
828: /* Find first boundary node */
829: for(corner = 0; corner < numCorners; corner++) {
830: MeshGetNodeFromElement(mesh, elem, corner, &minNode);
831: MeshGetNodeBoundary(mesh, minNode, &bd);
832: if (bd != 0) {
833: MeshGetNodeCoords(mesh, minNode, &nx, &ny, &nz);
834: minDist = PetscSqr(MeshPeriodicDiffX(mesh, nx - x)) + PetscSqr(MeshPeriodicDiffY(mesh, ny - y));
835: break;
836: }
837: }
838: if (corner == numCorners) SETERRQ1(PETSC_ERR_ARG_WRONG, "Element %d has no node on a boundary", elem);
839: /* Find closest node with field defined */
840: for(corner = 1; corner < numCorners; corner++) {
841: MeshGetNodeFromElement(mesh, elem, corner, &nearNode);
842: MeshGetNodeCoords(mesh, nearNode, &nx, &ny, &nz);
843: MeshGetNodeBoundary(mesh, nearNode, &bd);
844: dist = PetscSqr(nx - x) + PetscSqr(ny - y);
845: if ((bd != 0) && (dist < minDist)) {
846: minDist = dist;
847: minNode = nearNode;
848: }
849: }
850: PartitionLocalToGlobalNodeIndex(p, minNode, &globalMinNode);
851: } else {
852: minNode = -1;
853: globalMinNode = -1;
854: }
856: /* The minimum node might still be a ghost node */
857: MPI_Allreduce(&globalMinNode, &allMinNode, 1, MPI_INT, MPI_MAX, mesh->comm);
858: if ((allMinNode >= startNode) && (allMinNode < endNode)) {
859: *node = allMinNode - startNode;
860: } else {
861: *node = -1;
862: }
863: if (allMinNode < 0) PetscFunctionReturn(1);
864: return(0);
865: }
867: /*@
868: MeshGetNodeSupport - This function get the canonical element numbers of all elements
869: within the support of a given basis function. A call to MeshRestoreNodeSupport() must
870: be made before another call to this function.
872: Not collective
874: Input Parameters:
875: + mesh - The mesh
876: . node - The node containing the basis function
877: - elem - An element containing the node, -1 for a search
879: Output Parameters:
880: + degree - The degree of the node
881: - support - A list of the elements in the support
883: Note:
884: This function currently only returns the elements containing
885: any given node, so some basis functions will have a wider
886: support than this definition.
888: Level: intermediate
890: .keywords: mesh, support
891: .seealso: MeshRestoreNodeSupport()
892: @*/
893: int MeshGetNodeSupport(Mesh mesh, int node, int elem, int *degree, int **support)
894: {
901: if (mesh->supportTaken == PETSC_TRUE) {
902: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Missing call to MeshRestoreNodeSupport()");
903: }
904: mesh->supportTaken = PETSC_TRUE;
905: (*mesh->ops->getnodesupport)(mesh, node, elem, degree, support);
906: return(0);
907: }
909: /*@
910: MeshRestoreNodeSupport - This function returns the support array from MeshGetNodeSupport().
912: Not collective
914: Input Parameters:
915: + mesh - The mesh
916: . node - The node containing the basis function
917: . elem - An element containing the node, -1 for a search
918: . degree - The degree pf the node
919: - support - A list of the elements in the support
921: Level: intermediate
923: .keywords: mesh, support
924: .seealso: MeshGetNodeSupport()
925: @*/
926: int MeshRestoreNodeSupport(Mesh mesh, int node, int elem, int *degree, int **support)
927: {
930: if (*support != mesh->support) {
931: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid support argument");
932: }
933: mesh->supportTaken = PETSC_FALSE;
934: return(0);
935: }
937: /*@
938: MeshGetNodeOrdering - This function gets the AO which was used to reorder the mesh nodes
939: before partitioning. This is sometimes used to bandwdith reduce the mesh graph.
941: Not collective
943: Input Parameter:
944: . mesh - The mesh
946: Output Parameter:
947: . order - The node ordering, or PETSC_NULL if none exists
949: Level: intermediate
951: .keywords: Mesh, node, ordering, AO
952: .seealso: PartitionGetNodeOrdering()
953: @*/
954: int MeshGetNodeOrdering(Mesh mesh, AO *order)
955: {
959: *order = mesh->nodeOrdering;
960: return(0);
961: }
963: /*------------------------------------------- Mesh Element Query Functions -------------------------------------------*/
964: /*@
965: MeshGetElementNeighbor - This function retrieves the local element number of the element, or neighbor, opposing a
966: given corner.
968: Not collective
970: Input Parameters:
971: + mesh - The mesh
972: . elem - The element
973: - corner - The corner, or element node number
975: Output Parameter:
976: . neighbor - The neighbor opposing corner
978: Level: intermediate
980: .keywords: mesh, element, neighbor
981: .seealso: MeshGetNodeCoords()
982: @*/
983: int MeshGetElementNeighbor(Mesh mesh, int elem, int corner, int *neighbor)
984: {
990: (*mesh->ops->getelemneighbor)(mesh, elem, corner, neighbor);
991: return(0);
992: }
994: /*@
995: MeshLocatePoint - This function returns the element containing
996: the given point, or -1 if it is not contained in the local mesh.
998: Not collective
1000: Input Parameters:
1001: + mesh - The mesh
1002: - x,y,z - The node coordinates
1004: Output Parameter:
1005: . containingElem - An element containing the node, -1 for failure
1007: Level: beginner
1009: .keywords: mesh, point location
1010: .seealso MeshGetNearestNode(), MeshGetNodeSupport()
1011: @*/
1012: int MeshLocatePoint(Mesh mesh, double x, double y, double z, int *containingElem)
1013: {
1019: PetscLogEventBegin(MESH_LocatePoint, mesh, 0, 0, 0);
1020: (*mesh->ops->locatepoint)(mesh, x, y, z, containingElem);
1021: PetscLogEventEnd(MESH_LocatePoint, mesh, 0, 0, 0);
1022: return(0);
1023: }
1025: /*--------------------------------------------- Mesh Hole Query Functions --------------------------------------------*/
1026: /*@
1027: MeshSetHoleCoords - Sets the coordinates of holes in the mesh
1029: Collective on Mesh
1031: Input Parameter:
1032: + mesh - The mesh
1033: . num - The number of holes
1034: - coords - The coordinates
1036: Level: advanced
1038: .keywords: mesh, hole, coordinates
1039: .seealso: MeshGetBoundaryStart()
1040: @*/
1041: int MeshSetHoleCoords(Mesh mesh, int num, Vec coords)
1042: {
1043: int dim = mesh->dim;
1044: PetscScalar *array;
1045: int size, vecDim;
1046: int h, d;
1047: int ierr;
1052: if (num != mesh->numHoles) {
1053: mesh->numHoles = num;
1054: if (mesh->holes != PETSC_NULL) {
1055: PetscFree(mesh->holes);
1056: }
1057: if (mesh->numHoles > 0) {
1058: PetscMalloc(mesh->numHoles*dim * sizeof(double), &mesh->holes);
1059: } else {
1060: mesh->holes = PETSC_NULL;
1061: }
1062: }
1063: VecGetArray(coords, &array);
1064: VecGetLocalSize(coords, &size);
1065: if (size == mesh->numHoles*dim) {
1066: PetscMemcpy(mesh->holes, array, mesh->numHoles*dim * sizeof(double));
1067: } else if ((size > mesh->numHoles*dim) && (size%mesh->numHoles == 0)) {
1068: vecDim = (int) (size/mesh->numHoles);
1069: for(h = 0; h < mesh->numHoles; h++) {
1070: for(d = 0; d < dim; d++) {
1071: mesh->holes[h*dim+d] = PetscRealPart(array[h*vecDim+d]);
1072: }
1073: }
1074: } else {
1075: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid vector of hole coordinates");
1076: }
1077: VecRestoreArray(coords, &array);
1078: return(0);
1079: }
1081: /*------------------------------------------ Mesh Embedding Query Functions ------------------------------------------*/
1082: /*@
1083: MeshGetElementFromNode - This function retrieves a local element number from the local
1084: node number. Notice that this could be nondeterministic.
1086: Not collective
1088: Input Parameters:
1089: + mesh - The mesh
1090: - node - The number
1092: Output Parameter:
1093: . elem - The element
1095: Level: intermediate
1097: .keywords: mesh, node, element
1098: .seealso: MeshGetNodeFromElement()
1099: @*/
1100: int MeshGetElementFromNode(Mesh mesh, int node, int *elem)
1101: {
1102: Partition part;
1103: int numElements, numCorners;
1104: int e, c, n;
1105: int ierr;
1110: MeshGetPartition(mesh, &part);
1111: PartitionGetNumElements(part, &numElements);
1112: MeshGetNumCorners(mesh, &numCorners);
1113: for(e = 0; e < numElements; e++) {
1114: for(c = 0; c < numCorners; c++){
1115: MeshGetNodeFromElement(mesh, e, c, &n);
1116: if (n == node) {
1117: *elem = e;
1118: return(0);
1119: }
1120: }
1121: }
1122: PetscFunctionReturn(-1);
1123: }
1125: /*@
1126: MeshGetBdElementFromEdge - This function retrieves the local element number of a
1127: boundary element from the local edge number.
1129: Not collective
1131: Input Parameters:
1132: + mesh - The mesh
1133: - edge - The edge
1135: Output Parameter:
1136: . elem - The element along the boundary
1138: Level: intermediate
1140: .keywords: mesh, element, edge
1141: .seealso: MeshGetNodeFromElement(), MeshGetNodeCoords()
1142: @*/
1143: int MeshGetBdElementFromEdge(Mesh mesh, int edge, int *elem)
1144: {
1145: #if 0
1146: int neighbor;
1147: #endif
1148: int degree;
1149: int *support;
1150: int startNode, endNode;
1151: int startCorner, endCorner, numCorners;
1152: int sElem, corner, node;
1153: int ierr;
1157: /* This should be a default function, which is only used if the implementation
1158: does not have an optimized implementation
1159: */
1160: /* Search for element containing edge */
1161: MeshGetNumCorners(mesh, &numCorners);
1162: MeshGetNodeFromEdge(mesh, edge, 0, &startNode);
1163: MeshGetNodeFromEdge(mesh, edge, 1, &endNode);
1164: /* Loop over nodes on each element in the support of the node */
1165: MeshGetNodeSupport(mesh, startNode, 0, °ree, &support);
1166: for(sElem = 0; sElem < degree; sElem++) {
1167: for(corner = 0, startCorner = -1, endCorner = -1; corner < numCorners; corner++) {
1168: MeshGetNodeFromElement(mesh, support[sElem], corner, &node);
1169: if (node == startNode) {
1170: startCorner = corner;
1171: } else if (node == endNode) {
1172: endCorner = corner;
1173: }
1174: }
1175: if ((startCorner >= 0) && (endCorner >= 0)) break;
1176: }
1177: if (sElem == degree) SETERRQ1(PETSC_ERR_PLIB, "Edge %d not found in element list", edge);
1178: *elem = support[sElem];
1179: MeshRestoreNodeSupport(mesh, startNode, 0, °ree, &support);
1180: #if 0
1181: MeshGetNeighbor(mesh, *elem, ((startCorner+endCorner)*2)%3, &neighbor);
1182: if (neighbor != -1)
1183: SETERRQ2(PETSC_ERR_ARG_WRONG, "Edge is not on a boundary since elem %d has neighbor %d", *elem, neighbor);
1184: #endif
1185: return(0);
1186: }
1188: /*@
1189: MeshGetNodeFromElement - This function retrieves the local node number from the local
1190: element number and the corner.
1192: Not collective
1194: Input Parameters:
1195: + mesh - The mesh
1196: . elem - The element
1197: - corner - The corner, or element node number
1199: Output Parameter:
1200: . node - The local node number
1202: Level: intermediate
1204: .keywords: mesh, node
1205: .seealso: MeshGetNodeCoords()
1206: @*/
1207: int MeshGetNodeFromElement(Mesh mesh, int elem, int corner, int *node)
1208: {
1209: int numCorners = mesh->numCorners;
1210: int numOverlapElements = mesh->numFaces;
1216: if (mesh->part != PETSC_NULL) {
1217: ierr = PartitionGetNumOverlapElements(mesh->part, &numOverlapElements);
1218: }
1219: if ((elem < 0) || (elem >= numOverlapElements)) {
1220: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node %d should be in [0,%d)", elem, numOverlapElements);
1221: }
1222: if ((corner < 0) || (corner >= numCorners)) {
1223: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid corner %d should be in [0,%d)", elem, numCorners);
1224: }
1225: (*mesh->ops->getnodefromelement)(mesh, elem, corner, node);
1226: return(0);
1227: }
1229: /*@
1230: MeshGetNodeFromEdge - This function retrieves the local node number from the local
1231: edge number and the corner.
1233: Not collective
1235: Input Parameters:
1236: + mesh - The mesh
1237: . edge - The edge
1238: - corner - The corner, or edge node number (0 or 1)
1240: Output Parameter:
1241: . node - The local node number
1243: Level: intermediate
1245: .keywords: mesh, node, edge
1246: .seealso: MeshGetNodeCoords()
1247: @*/
1248: int MeshGetNodeFromEdge(Mesh mesh, int edge, int corner, int *node)
1249: {
1250: int numOverlapEdges = mesh->numEdges;
1256: if ((corner < 0) || (corner >= 2)) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Invalid corner %d should be in [0,2)", edge);
1257: if (mesh->part != PETSC_NULL) {
1258: ierr = PartitionGetNumOverlapEdges(mesh->part, &numOverlapEdges);
1259: }
1260: if ((edge < 0) || (edge >= numOverlapEdges)) {
1261: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid edge %d should be in [0,%d)", edge, numOverlapEdges);
1262: }
1263: (*mesh->ops->getnodefromedge)(mesh, edge, corner, node);
1264: return(0);
1265: }
1267: /*@
1268: MeshGetMidnodeFromEdge - This function retrieves the local node number of the midnode from the local edge number.
1270: Not collective
1272: Input Parameters:
1273: + mesh - The mesh
1274: - edge - The edge
1276: Output Parameter:
1277: . midnode - The local node number of the midnode
1279: Level: intermediate
1281: .keywords: mesh, midnode, edge
1282: .seealso: MeshGetNodeFromElement(), MeshGetNodeCoords()
1283: @*/
1284: int MeshGetMidnodeFromEdge(Mesh mesh, int edge, int *midnode)
1285: {
1286: int degree;
1287: int *support;
1288: int startNode, endNode;
1289: int numCorners;
1290: int startCorner = -1;
1291: int endCorner = -1;
1292: int elem, corner, node;
1293: int ierr;
1297: /* This should be a default function, which is only used if the implementation
1298: does not have an optimized implementation
1299: */
1300: MeshGetNumCorners(mesh, &numCorners);
1301: if (numCorners == 3) return(0);
1302: if (numCorners != 6) SETERRQ1(PETSC_ERR_ARG_WRONG, "Cannot find midnode for %d corners", numCorners);
1304: /* Search for midnode on edge */
1305: MeshGetNodeFromEdge(mesh, edge, 0, &startNode);
1306: MeshGetNodeFromEdge(mesh, edge, 1, &endNode);
1307: /* Loop over nodes on each element in the support of the node */
1308: MeshGetNodeSupport(mesh, startNode, 0, °ree, &support);
1309: for(elem = 0; elem < degree; elem++) {
1310: for(corner = 0, startCorner = -1, endCorner = -1; corner < numCorners; corner++) {
1311: MeshGetNodeFromElement(mesh, support[elem], corner, &node);
1312: if (node == startNode) {
1313: startCorner = corner;
1314: } else if (node == endNode) {
1315: endCorner = corner;
1316: }
1317: }
1318: if ((startCorner >= 0) && (endCorner >= 0)) break;
1319: }
1320: if (elem == degree) {
1321: if (startCorner < 0) SETERRQ1(PETSC_ERR_PLIB, "Invalid support of node %d did not contain node", startNode);
1322: if (endCorner < 0) {
1323: for(elem = 0; elem < mesh->numFaces; elem++) {
1324: for(corner = 0, startCorner = -1, endCorner = -1; corner < numCorners; corner++) {
1325: MeshGetNodeFromElement(mesh, elem, corner, &node);
1326: if (node == startNode) {
1327: startCorner = corner;
1328: } else if (node == endNode) {
1329: endCorner = corner;
1330: }
1331: }
1332: if ((startCorner >= 0) && (endCorner >= 0)) break;
1333: }
1334: if (elem == mesh->numFaces) {
1335: SETERRQ3(PETSC_ERR_PLIB, "Edge %d (%d,%d) not found in element list", edge, startNode, endNode);
1336: } else {
1337: SETERRQ5(PETSC_ERR_PLIB, "Edge %d (%d,%d) found in element %d not in support list of node %d",
1338: edge, startNode, endNode, elem, startNode);
1339: }
1340: }
1341: }
1343: MeshGetNodeFromElement(mesh, support[elem], ((startCorner + endCorner)*2)%3 + 3, midnode);
1344: MeshRestoreNodeSupport(mesh, startNode, 0, °ree, &support);
1345: return(0);
1346: }