Actual source code: tri2dMovement.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: tri2dMovement.c,v 1.4 2000/01/10 03:54:23 knepley Exp $";
3: #endif
5: /* Mesh movement for 2d triangular grids */
6: #include "src/gvec/meshMovementImpl.h" /*I "meshMovement.h" I*/
7: #include "src/mesh/impls/triangular/2d/2dimpl.h" /*I "mesh.h" I*/
8: #include "src/gvec/gvecimpl.h" /* For MeshMoveMesh() */
9: #include "tri2dMovement.h"
11: int MeshMoverSetup_Triangular_2D(MeshMover mover) {
12: Mesh mesh = mover->mesh;
13: Mesh movingMesh = mesh;
14: int dim = mesh->dim;
15: int numCorners = mesh->numCorners;
16: DiscretizationType dtype;
17: PetscTruth opt;
18: int ierr;
21: /* Determine order of interpolation */
22: if (numCorners == 3) {
23: dtype = DISCRETIZATION_TRIANGULAR_2D_LINEAR;
24: } else if (numCorners == 6) {
25: dtype = DISCRETIZATION_TRIANGULAR_2D_QUADRATIC;
26: ierr = PetscOptionsHasName(mesh->prefix, "-mesh_move_linear", &opt);
27: if (opt == PETSC_TRUE) {
28: dtype = DISCRETIZATION_TRIANGULAR_2D_LINEAR;
29: /* Must play games here so that the coarse mesh does not get a MeshMover as well */
30: PetscObjectReference((PetscObject) mover);
31: MeshSetMover(mesh, PETSC_NULL);
32: MeshCoarsen(mesh, PETSC_NULL, &movingMesh);
33: MeshSetMover(mesh, mover);
34: }
35: } else {
36: SETERRQ1(PETSC_ERR_SUP, "Number of nodes per element %d not supported", numCorners);
37: }
39: /* Setup mesh velocity grid */
40: GridCreate(movingMesh, &mover->meshVelocityGrid);
41: GridAddField(mover->meshVelocityGrid, "Velocity", dtype, dim, PETSC_NULL);
42: GridAppendOptionsPrefix(mover->meshVelocityGrid, "mesh_vel_");
43: GridSetFromOptions(mover->meshVelocityGrid);
44: GridSetActiveField(mover->meshVelocityGrid, 0);
45: PetscLogObjectParent(movingMesh, mover->meshVelocityGrid);
46: /* ALE operators need boundary values */
47: mover->meshVelocityGrid->reduceElement = PETSC_TRUE;
48: switch(mover->meshVelocityMethod) {
49: case MESH_LAPLACIAN:
50: GridSetMatOperator(mover->meshVelocityGrid, 0, 0, LAPLACIAN, 1.0, PETSC_FALSE, PETSC_NULL);
51: break;
52: case MESH_WEIGHTED_LAPLACIAN:
53: GridSetMatOperator(mover->meshVelocityGrid, 0, 0, WEIGHTED_LAP, 1.0, PETSC_FALSE, PETSC_NULL);
54: break;
55: default:
56: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid velocity calculation method %d", mover->meshVelocityMethod);
57: }
59: /* Setup mesh acceleration grid */
60: GridCreate(movingMesh, &mover->meshAccelerationGrid);
61: GridAddField(mover->meshAccelerationGrid, "Acceleration", dtype, dim, PETSC_NULL);
62: GridAppendOptionsPrefix(mover->meshAccelerationGrid, "mesh_accel_");
63: GridSetFromOptions(mover->meshAccelerationGrid);
64: GridSetActiveField(mover->meshAccelerationGrid, 0);
65: PetscLogObjectParent(mesh, mover->meshAccelerationGrid);
66: switch(mover->meshAccelerationMethod) {
67: case MESH_LAPLACIAN:
68: GridSetMatOperator(mover->meshAccelerationGrid, 0, 0, LAPLACIAN, 1.0, PETSC_FALSE, PETSC_NULL);
69: break;
70: case MESH_WEIGHTED_LAPLACIAN:
71: GridSetMatOperator(mover->meshVelocityGrid, 0, 0, WEIGHTED_LAP, 1.0, PETSC_FALSE, PETSC_NULL);
72: break;
73: default:
74: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid acceleration calculation method%d", mover->meshAccelerationMethod);
75: }
76: /* This will cause the coarse mesh to be destroyed when no longer needed */
77: if (mesh != movingMesh) {
78: MeshDestroy(movingMesh);
79: }
80: return(0);
81: }
83: int MeshMoverDestroy_Triangular_2D(MeshMover mover) {
85: return(0);
86: }
88: int MeshMapVertices_Triangular_2D(GVec coarseVec, Mesh mesh, Vec vec) {
89: Partition p = mesh->part;
90: int dim = mesh->dim;
91: Grid coarseGrid;
92: Mesh coarseMesh;
93: Partition coarseP;
94: AO coarseMap;
95: PetscTruth reduceSystem;
96: int numCoarseNodes;
97: int *classes;
98: int firstVar;
99: int *offsets, *localOffsets;
100: int reduceFirstVar = 0;
101: int *reduceOffsets = PETSC_NULL;
102: int *reduceLocalOffsets = PETSC_NULL;
103: int **reduceFieldClasses = PETSC_NULL;
104: PetscScalar *vals, *coarseVals, *reduceVals;
105: int node, cNode, nclass, coarseVar, j;
106: int ierr;
109: GVecGetGrid(coarseVec, &coarseGrid);
110: GridGetMesh(coarseGrid, &coarseMesh);
111: VecGetArray(vec, &vals);
112: VecGetArray(coarseVec, &coarseVals);
113: coarseP = coarseMesh->part;
114: coarseMap = coarseMesh->coarseMap;
115: reduceSystem = coarseGrid->reduceSystem;
116: numCoarseNodes = coarseMesh->numNodes;
117: reduceFieldClasses = coarseGrid->cm->reduceFieldClasses;
118: classes = coarseGrid->cm->classes;
119: firstVar = coarseGrid->order->firstVar[p->rank];
120: offsets = coarseGrid->order->offsets;
121: localOffsets = coarseGrid->order->localOffsets;
122: if (reduceSystem == PETSC_TRUE) {
123: reduceOffsets = coarseGrid->reduceOrder->offsets;
124: reduceLocalOffsets = coarseGrid->reduceOrder->localOffsets;
125: reduceFirstVar = coarseGrid->reduceOrder->firstVar[p->rank];
126: VecGetArray(coarseGrid->bdReduceVecCur, &reduceVals);
127: }
128: for(cNode = 0; cNode < numCoarseNodes; cNode++) {
129: /* Get corresponding node */
130: PartitionLocalToGlobalNodeIndex(coarseP, cNode, &node);
131: AOPetscToApplication(coarseMap, 1, &node);
132: PartitionGlobalToLocalNodeIndex(p, node, &node);
133: /* Calculate coarse offset */
134: nclass = classes[cNode];
135: if ((reduceSystem == PETSC_TRUE) && (reduceFieldClasses[0][nclass])) {
136: if (cNode >= numCoarseNodes) {
137: coarseVar = reduceLocalOffsets[cNode-numCoarseNodes];
138: } else {
139: coarseVar = reduceOffsets[cNode] - reduceFirstVar;
140: }
141: /* Set node value */
142: for(j = 0; j < dim; j++) {
143: vals[node*dim+j] = reduceVals[coarseVar+j];
144: }
145: } else {
146: if (cNode >= numCoarseNodes) {
147: coarseVar = localOffsets[cNode-numCoarseNodes];
148: } else {
149: coarseVar = offsets[cNode] - firstVar;
150: }
151: /* Set node value */
152: for(j = 0; j < dim; j++) {
153: vals[node*dim+j] = coarseVals[coarseVar+j];
154: }
155: }
156: }
157: VecRestoreArray(vec, &vals);
158: VecRestoreArray(coarseVec, &coarseVals);
159: if (reduceSystem == PETSC_TRUE) {
160: VecRestoreArray(coarseGrid->bdReduceVecCur, &reduceVals);
161: }
162: return(0);
163: }
165: int MeshInterpMidnodes_Triangular_2D(Mesh mesh, Vec vec) {
166: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
167: int dim = mesh->dim;
168: int numCorners = mesh->numCorners;
169: int numElements = mesh->numFaces;
170: int *elements = tri->faces;
171: PetscScalar *vals;
172: int elem, corner, node, nNode1, nNode2, j;
173: int ierr;
176: VecGetArray(vec, &vals);
177: for(elem = 0; elem < numElements; elem++) {
178: for(corner = 3; corner < numCorners; corner++) {
179: node = elements[elem*numCorners+corner];
180: nNode1 = elements[elem*numCorners+((corner+1)%3)];
181: nNode2 = elements[elem*numCorners+((corner+2)%3)];
182: for(j = 0; j < dim; j++) {
183: vals[node*dim+j] = 0.5*(vals[nNode1*dim+j] + vals[nNode2*dim+j]);
184: }
185: }
186: }
187: VecRestoreArray(vec, &vals);
188: PetscLogFlops(2*dim*(numCorners-3)*numElements);
189: return(0);
190: }
192: int MeshUpdateNodeValues_Triangular_2D(MeshMover mover, GVec sol, Vec vec, Vec ghostVec) {
193: Mesh mesh = mover->mesh;
194: Partition p = mesh->part;
195: int numCorners = mesh->numCorners;
196: Grid grid;
197: int numOverlapNodes, numFuncs;
198: PetscTruth match;
199: int ierr;
202: /* Check grid */
203: GVecGetGrid(sol, &grid);
204: PartitionGetNumOverlapNodes(p, &numOverlapNodes);
205: DiscretizationGetNumFunctions(grid->fields[0].disc, &numFuncs);
206: if (grid->reduceSystem == PETSC_TRUE) {
207: if (numFuncs == numCorners) {
208: if (grid->order->numOverlapVars + grid->reduceOrder->numOverlapVars != numOverlapNodes*2) {
209: SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid number of mesh velocity unknowns %d should be %d",
210: grid->order->numOverlapVars + grid->reduceOrder->numOverlapVars, numOverlapNodes*2);
211: }
212: }
213: } else {
214: if (numFuncs == numCorners) {
215: if (grid->order->numOverlapVars != numOverlapNodes*2) {
216: SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid number of mesh velocity unknowns %d should be %d",
217: grid->order->numOverlapVars, numOverlapNodes*2);
218: }
219: }
220: }
221: /* Map coarse nodes to fine vertices */
222: PetscTypeCompare((PetscObject) grid->fields[0].disc, DISCRETIZATION_TRIANGULAR_2D_LINEAR, &match);
223: if ((numCorners == 6) && (match == PETSC_TRUE)) {
224: MeshMapVertices_Triangular_2D(sol, mesh, vec);
225: } else {
226: VecCopy(sol, vec);
227: }
228: /* Copy to ghost vector */
229: VecScatterBegin(vec, ghostVec, INSERT_VALUES, SCATTER_FORWARD, mover->meshVelocityScatter);
230: VecScatterEnd(vec, ghostVec, INSERT_VALUES, SCATTER_FORWARD, mover->meshVelocityScatter);
231: /* Interpolate velocities of midnodes */
232: MeshInterpMidnodes_Triangular_2D(mesh, ghostVec);
234: /* REMOVE THIS AFTER YOU REWRITE THE ALE STUFF */
235: if (sol == mover->meshVelocitySol) {
236: GridGlobalToLocal(mover->meshVelocityGrid, INSERT_VALUES, sol);
237: } else if (sol == mover->meshAccelerationSol) {
238: GridGlobalToLocal(mover->meshAccelerationGrid, INSERT_VALUES, sol);
239: } else {
240: SETERRQ(PETSC_ERR_SUP, "Unrecognized solution vector");
241: }
242: return(0);
243: }
245: int MeshUpdateCoords_Triangular_2D(Mesh mesh, PetscScalar *v, PetscScalar *a, double t) {
246: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
247: Partition p = mesh->part;
248: int dim = mesh->dim;
249: double *nodes = tri->nodes;
250: int numOverlapNodes;
251: int node;
252: int ierr;
255: PartitionGetNumOverlapNodes(p, &numOverlapNodes);
256: /* x' = x + v t + 0.5 a t^2 */
257: if (mesh->isPeriodic == PETSC_TRUE) {
258: for(node = 0; node < numOverlapNodes; node++) {
259: nodes[node*2] = MeshPeriodicX(mesh, nodes[node*2] + PetscRealPart(v[node*2])*t + 0.5*PetscRealPart(a[node*2])*t*t);
260: nodes[node*2+1] = MeshPeriodicY(mesh, nodes[node*2+1] + PetscRealPart(v[node*2+1])*t + 0.5*PetscRealPart(a[node*2+1])*t*t);
261: #ifdef DEBUG_BD_MOVEMENT
262: if (tri->markers[node] > 0)
263: PetscPrintf(PETSC_COMM_SELF, "node %d: u = %lf v = %lfn", node, v[node*2], v[node*2+1]);
264: #endif
265: }
266: } else {
267: for(node = 0; node < numOverlapNodes*dim; node++) {
268: nodes[node] += PetscRealPart(v[node])*t + 0.5*PetscRealPart(a[node])*t*t;
269: #ifdef DEBUG_BD_MOVEMENT
270: if ((node%2 == 0) && (tri->markers[node/2] > 0))
271: PetscPrintf(PETSC_COMM_SELF, "node %d: u = %lf v = %lfn", node/2, v[node], v[node+1]);
272: #endif
273: }
274: }
275: PetscLogFlops(6*numOverlapNodes*dim);
276: return(0);
277: }
279: int MeshSyncCoords_Triangular_2D(Mesh mesh, Mesh coarseMesh) {
280: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
281: Mesh_Triangular *coarseTri = (Mesh_Triangular *) coarseMesh->data;
282: Partition p = mesh->part;
283: Partition coarseP = coarseMesh->part;
284: AO coarseMap = coarseMesh->coarseMap;
285: int dim = mesh->dim;
286: int numCoarseNodes = coarseMesh->numNodes;
287: double *nodes = tri->nodes;
288: double *coarseNodes = coarseTri->nodes;
289: int node, cNode, j;
290: int ierr;
293: for(cNode = 0; cNode < numCoarseNodes; cNode++) {
294: PartitionLocalToGlobalNodeIndex(coarseP, cNode, &node);
295: AOPetscToApplication(coarseMap, 1, &node);
296: PartitionGlobalToLocalNodeIndex(p, node, &node);
297: for(j = 0; j < dim; j++) {
298: coarseNodes[cNode*dim+j] = nodes[node*dim+j];
299: }
300: }
301: return(0);
302: }
304: int MeshMoveMesh_Triangular_2D(MeshMover mover, double t) {
305: Mesh mesh = mover->mesh;
306: Mesh velMesh = mover->meshVelocityGrid->mesh;
307: PetscScalar *v, *a;
308: PetscDraw draw;
309: int ierr;
312: VecGetArray(mover->nodeVelocitiesGhost, &v);
313: VecGetArray(mover->nodeAccelerationsGhost, &a);
315: /* Move mesh nodes */
316: MeshUpdateCoords_Triangular_2D(mesh, v, a, t);
317: /* Move mesh velocity and acceleration mesh nodes */
318: MeshSyncCoords_Triangular_2D(mesh, velMesh);
320: VecRestoreArray(mover->nodeVelocitiesGhost, &v);
321: VecRestoreArray(mover->nodeAccelerationsGhost, &a);
322:
323: /* Update bounding rectangle */
324: MeshUpdateBoundingBox(mesh);
326: /* View mesh */
327: if (mover->movingMeshViewer) {
328: PetscViewerDrawGetDraw(mover->movingMeshViewer, 0, &draw);
329: PetscDrawClear(draw);
330: PetscDrawSetTitle(draw, mover->movingMeshViewerCaption);
331: MeshView(mesh, mover->movingMeshViewer);
332: PetscViewerFlush(mover->movingMeshViewer);
333: PetscDrawPause(draw);
334: }
335: return(0);
336: }
338: /*
339: MeshReform_Triangular_2D - Reforms a two dimensional unstructured mesh
340: using the original boundary.
342: Input Parameters:
343: . mesh - The mesh begin reformed
344: . refine - Flag indicating whether to refine or recreate the mesh
345: . area - A function which gives an area constraint when evaluated inside an element
346: . newBd - Flag to determine whether the boundary should be generated (PETSC_TRUE) or taken from storage
348: Output Parameter:
349: . newMesh - The reformed mesh
351: .keywords unstructured grid
352: .seealso GridReform(), MeshRefine(), MeshSetBoundary()
353: */
354: int MeshReform_Triangular_2D(Mesh mesh, PetscTruth refine, PointFunction area, PetscTruth newBd, Mesh *newMesh) {
355: Mesh m;
356: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
357: Partition p = mesh->part;
358: Partition_Triangular_2D *q = (Partition_Triangular_2D *) p->data;
359: MeshBoundary2D bdCtx;
360: int numEdges = mesh->numBdEdges;
361: int origLocNodes = mesh->numBdNodes;
362: int numLocNodes;
363: int *firstBdEdge;
364: int *firstBdNode;
365: int *firstHole;
366: int *numRecvNodes;
367: int *numRecvMarkers;
368: int *numRecvSegments;
369: int *numRecvSegMarkers;
370: int *numRecvHoles;
371: int *nodeOffsets;
372: int *markerOffsets;
373: int *segmentOffsets;
374: int *segMarkerOffsets;
375: int *holeOffsets;
376: double *locNodes;
377: int *locMarkers;
378: int *locSegments;
379: int *locSegMarkers;
380: double *nodes;
381: int *markers;
382: int *segments;
383: int *segMarkers;
384: double *holes;
385: int *globalNodes;
386: int *aOrder;
387: int rank, numProcs;
388: int numLocSegments, numSegments, numNodes, numHoles;
389: int proc, bd, gNode, lNode, bdNode, edge, bdEdge, edgeEnd;
390: int ierr;
393: MPI_Comm_size(mesh->comm, &numProcs);
394: MPI_Comm_rank(mesh->comm, &rank);
395: if ((newBd == PETSC_TRUE) || (mesh->bdCtxNew == PETSC_NULL)) {
396: MPI_Reduce(&mesh->numHoles, &numHoles, 1, MPI_INT, MPI_SUM, 0, mesh->comm);
397: PetscMalloc((numProcs+1) * sizeof(int), &firstBdEdge);
398: PetscMalloc((numProcs+1) * sizeof(int), &firstBdNode);
399: PetscMalloc((numProcs+1) * sizeof(int), &firstHole);
400: PetscMalloc(numProcs * sizeof(int), &numRecvNodes);
401: PetscMalloc(numProcs * sizeof(int), &numRecvMarkers);
402: PetscMalloc(numProcs * sizeof(int), &numRecvSegments);
403: PetscMalloc(numProcs * sizeof(int), &numRecvSegMarkers);
404: PetscMalloc(numProcs * sizeof(int), &numRecvHoles);
405: PetscMalloc(numProcs * sizeof(int), &nodeOffsets);
406: PetscMalloc(numProcs * sizeof(int), &markerOffsets);
407: PetscMalloc(numProcs * sizeof(int), &segmentOffsets);
408: PetscMalloc(numProcs * sizeof(int), &segMarkerOffsets);
409: PetscMalloc(numProcs * sizeof(int), &holeOffsets);
410: PetscMalloc(origLocNodes*2 * sizeof(double), &locNodes);
411: PetscMalloc(origLocNodes * sizeof(int), &locMarkers);
412: PetscMalloc(numEdges*2 * sizeof(int), &locSegments);
413: PetscMalloc(numEdges * sizeof(int), &locSegMarkers);
414: PetscMalloc(origLocNodes * sizeof(int), &globalNodes);
416: /* Make all node numbers global */
417: numLocNodes = 0;
418: numLocSegments = 0;
419: for(bd = 0; bd < mesh->numBd; bd++) {
420: /* Get boundary nodes */
421: for(bdNode = tri->bdBegin[bd]; bdNode < tri->bdBegin[bd+1]; bdNode++) {
422: gNode = tri->bdNodes[bdNode];
423: lNode = gNode - q->firstNode[p->rank];
424: if ((lNode >= 0) && (lNode < mesh->numNodes)) {
425: if (tri->degrees[lNode] > 0) {
426: locNodes[numLocNodes*2] = tri->nodes[lNode*2];
427: locNodes[numLocNodes*2+1] = tri->nodes[lNode*2+1];
428: locMarkers[numLocNodes] = tri->bdMarkers[bd];
429: globalNodes[numLocNodes] = gNode;
430: numLocNodes++;
431: }
432: }
433: }
435: /* Get boundary segments */
436: for(bdEdge = tri->bdEdgeBegin[bd]; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++) {
437: edge = tri->bdEdges[bdEdge] - q->firstEdge[rank];
438: if ((edge >= 0) && (edge < mesh->numEdges)) {
439: for(edgeEnd = 0; edgeEnd < 2; edgeEnd++) {
440: lNode = tri->edges[edge*2+edgeEnd];
441: ierr = PartitionLocalToGlobalNodeIndex(p, lNode, &gNode);
442: locSegments[numLocSegments*2+edgeEnd] = gNode;
443: }
444: locSegMarkers[numLocSegments] = tri->bdMarkers[bd];
445: numLocSegments++;
446: }
447: }
448: }
450: /* Assemble offsets */
451: MPI_Allgather(&numLocSegments, 1, MPI_INT, &firstBdEdge[1], 1, MPI_INT, p->comm);
452: firstBdEdge[0] = 0;
453: for(proc = 1; proc <= numProcs; proc++)
454: firstBdEdge[proc] = firstBdEdge[proc] + firstBdEdge[proc-1];
455: numSegments = firstBdEdge[numProcs];
456: if ((rank == 0) && (numSegments != numEdges)) SETERRQ(PETSC_ERR_PLIB, "Invalid number of boundary segments");
457: MPI_Allgather(&numLocNodes, 1, MPI_INT, &firstBdNode[1], 1, MPI_INT, p->comm);
458: firstBdNode[0] = 0;
459: for(proc = 1; proc <= numProcs; proc++)
460: firstBdNode[proc] = firstBdNode[proc] + firstBdNode[proc-1];
461: numNodes = firstBdNode[numProcs];
462: MPI_Allgather(&mesh->numHoles, 1, MPI_INT, &firstHole[1], 1, MPI_INT, p->comm);
463: firstHole[0] = 0;
464: for(proc = 1; proc <= numProcs; proc++)
465: firstHole[proc] = firstHole[proc] + firstHole[proc-1];
466: if ((rank == 0) && (numHoles != firstHole[numProcs])) SETERRQ(PETSC_ERR_PLIB, "Invalid number of holes");
468: for(proc = 0; proc < numProcs; proc++) {
469: numRecvNodes[proc] = (firstBdNode[proc+1] - firstBdNode[proc])*2;
470: numRecvMarkers[proc] = (firstBdNode[proc+1] - firstBdNode[proc]);
471: numRecvSegments[proc] = (firstBdEdge[proc+1] - firstBdEdge[proc])*2;
472: numRecvSegMarkers[proc] = (firstBdEdge[proc+1] - firstBdEdge[proc]);
473: numRecvHoles[proc] = (firstHole[proc+1] - firstHole[proc])*2;
474: }
475: nodeOffsets[0] = 0;
476: markerOffsets[0] = 0;
477: segmentOffsets[0] = 0;
478: segMarkerOffsets[0] = 0;
479: holeOffsets[0] = 0;
480: for(proc = 1; proc < numProcs; proc++) {
481: nodeOffsets[proc] = numRecvNodes[proc-1] + nodeOffsets[proc-1];
482: markerOffsets[proc] = numRecvMarkers[proc-1] + markerOffsets[proc-1];
483: segmentOffsets[proc] = numRecvSegments[proc-1] + segmentOffsets[proc-1];
484: segMarkerOffsets[proc] = numRecvSegMarkers[proc-1] + segMarkerOffsets[proc-1];
485: holeOffsets[proc] = numRecvHoles[proc-1] + holeOffsets[proc-1];
486: }
488: if (rank == 0) {
489: PetscMalloc(numNodes*2 * sizeof(double), &nodes);
490: PetscMalloc(numNodes * sizeof(int), &markers);
491: PetscMalloc(numEdges*2 * sizeof(int), &segments);
492: PetscMalloc(numEdges * sizeof(int), &segMarkers);
493: if (numHoles) {
494: ierr = PetscMalloc(numHoles*2 * sizeof(double), &holes);
495: } else {
496: holes = PETSC_NULL;
497: }
498: }
500: /* Make new numbering for boundary nodes */
501: PetscMalloc(numNodes * sizeof(int), &aOrder);
502: MPI_Allgatherv(globalNodes, numLocNodes, MPI_INT, aOrder, numRecvMarkers, markerOffsets, MPI_INT, p->comm);
503:
505: /* Renumber segments */
506: for(edge = 0; edge < numLocSegments*2; edge++) {
507: gNode = locSegments[edge];
508: for(bdNode = 0; bdNode < numNodes; bdNode++)
509: if (aOrder[bdNode] == gNode)
510: break;
511: if (bdNode == numNodes) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary node");
512: locSegments[edge] = bdNode;
513: }
514: PetscFree(aOrder);
516: MPI_Gatherv(locNodes, numLocNodes*2, MPI_DOUBLE, nodes, numRecvNodes, nodeOffsets, MPI_DOUBLE, 0, mesh->comm);
517:
518: MPI_Gatherv(locMarkers, numLocNodes, MPI_INT, markers, numRecvMarkers, markerOffsets, MPI_INT, 0, mesh->comm);
519:
520: MPI_Gatherv(locSegments, numLocSegments*2, MPI_INT, segments, numRecvSegments, segmentOffsets, MPI_INT, 0, mesh->comm);
521:
522: MPI_Gatherv(locSegMarkers, numLocSegments, MPI_INT, segMarkers, numRecvSegMarkers, segMarkerOffsets, MPI_INT, 0, mesh->comm);
523:
524: MPI_Gatherv(mesh->holes, mesh->numHoles*2, MPI_DOUBLE, holes, numRecvHoles, holeOffsets, MPI_DOUBLE, 0, mesh->comm);
525:
527: /* Cleanup */
528: PetscFree(firstBdEdge);
529: PetscFree(firstBdNode);
530: PetscFree(firstHole);
531: PetscFree(numRecvNodes);
532: PetscFree(numRecvMarkers);
533: PetscFree(numRecvSegments);
534: PetscFree(numRecvSegMarkers);
535: PetscFree(numRecvHoles);
536: PetscFree(nodeOffsets);
537: PetscFree(markerOffsets);
538: PetscFree(segmentOffsets);
539: PetscFree(segMarkerOffsets);
540: PetscFree(holeOffsets);
541: PetscFree(locNodes);
542: PetscFree(locMarkers);
543: PetscFree(locSegments);
544: PetscFree(locSegMarkers);
545: PetscFree(globalNodes);
547: bdCtx.numVertices = numNodes;
548: bdCtx.vertices = nodes;
549: bdCtx.markers = markers;
550: bdCtx.numSegments = numEdges;
551: bdCtx.segments = segments;
552: bdCtx.segMarkers = segMarkers;
553: bdCtx.numHoles = numHoles;
554: bdCtx.holes = holes;
555: } else {
556: bdCtx.numVertices = mesh->bdCtxNew->numVertices;
557: bdCtx.vertices = mesh->bdCtxNew->vertices;
558: bdCtx.markers = mesh->bdCtxNew->markers;
559: bdCtx.numSegments = mesh->bdCtxNew->numSegments;
560: bdCtx.segments = mesh->bdCtxNew->segments;
561: bdCtx.segMarkers = mesh->bdCtxNew->segMarkers;
562: bdCtx.numHoles = mesh->bdCtxNew->numHoles;
563: bdCtx.holes = mesh->bdCtxNew->holes;
564: }
565: bdCtx.numBd = mesh->numBd;
567: MeshCreate(mesh->comm, &m);
568: MeshSetDimension(m, 2);
569: MeshSetPeriodicDimension(m, 0, mesh->isPeriodicDim[0]);
570: MeshSetPeriodicDimension(m, 1, mesh->isPeriodicDim[1]);
571: MeshSetNumCorners(m, mesh->numCorners);
572: MeshSetBoundary(m, &bdCtx);
573: MeshSetType(m, mesh->type_name);
574: PetscObjectSetName((PetscObject) m, "Mesh");
575: *newMesh = m;
577: if (rank == 0) {
578: if ((newBd == PETSC_TRUE) || (mesh->bdCtxNew == PETSC_NULL)) {
579: PetscFree(nodes);
580: PetscFree(markers);
581: PetscFree(segments);
582: PetscFree(segMarkers);
583: if (numHoles) {
584: PetscFree(holes);
585: }
586: } else {
587: #if 0
588: MeshBoundaryDestroy(mesh->bdCtxNew);
589: #endif
590: mesh->bdCtxNew = PETSC_NULL;
591: }
592: }
594: #ifdef PETSC_USE_BOPT_g
595: #endif
596: return(0);
597: }
599: int MeshDebugBoundary_Triangular_2D_Private(Mesh mesh, GVec nodeVelocities) {
600: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
601: int numNodes = mesh->numNodes;
602: double *nodes = tri->nodes;
603: int *markers = tri->markers;
604: int rank = mesh->part->rank;
605: PetscScalar *vel;
606: int node;
607: int ierr;
610: VecGetArray(nodeVelocities, &vel);
611: for(node = 0; node < numNodes; node++) {
612: if (markers[node] > 0) {
613: PetscPrintf(PETSC_COMM_SELF, "[%d]node %d at (%g,%g) vel: (%g,%g)n",
614: rank, node, nodes[node*2], nodes[node*2+1], PetscRealPart(vel[node*2]), PetscRealPart(vel[node*2+1]));
615: }
616: }
617: VecRestoreArray(nodeVelocities, &vel);
618: return(0);
619: }
621: int MeshDebugBoundary_Triangular_2D(MeshMover mover) {
622: Mesh mesh = mover->mesh;
623: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
624: Grid velGrid = mover->meshVelocityGrid;
625: Grid accGrid = mover->meshAccelerationGrid;
626: int numBd = mesh->numBd;
627: int *bdMarkers = tri->bdMarkers;
628: int numNodes = mesh->numNodes;
629: int rank = mesh->part->rank;
630: int **reduceFieldClasses = velGrid->cm->reduceFieldClasses;
631: int *classes = velGrid->cm->classes;
632: int firstVar = velGrid->order->firstVar[rank];
633: int *offsets = velGrid->order->offsets;
634: int *localOffsets = velGrid->order->localOffsets;
635: int reduceFirstVar = 0;
636: int *reduceOffsets = PETSC_NULL;
637: int *reduceLocalOffsets = PETSC_NULL;
638: PetscScalar *v, *reduceV, *reduceA;
639: int b, bd, node, nclass, var;
640: int ierr;
643: VecGetArray(velGrid->ghostVec, &v);
644: #ifdef DEBUG_ACC
645: VecGetArray(accGrid->ghostVec, &a);
646: #endif
647: if (velGrid->reduceSystem == PETSC_TRUE) {
648: reduceOffsets = velGrid->reduceOrder->offsets;
649: reduceLocalOffsets = velGrid->reduceOrder->localOffsets;
650: reduceFirstVar = velGrid->reduceOrder->firstVar[rank];
651: VecGetArray(velGrid->bdReduceVecCur, &reduceV);
652: VecGetArray(accGrid->bdReduceVecCur, &reduceA);
653: }
654: for(b = 0; b < numBd; b++) {
655: bd = bdMarkers[b];
656: if (bd < 0) continue;
658: MeshGetBoundaryStart(mesh, bd, PETSC_TRUE, &node);
659: while (node >= 0) {
660: /* Print the value of velocity on boundary node */
661: nclass = classes[node];
662: if ((reduceFieldClasses != PETSC_NULL) && (reduceFieldClasses[0][nclass])) {
663: if (node >= numNodes) {
664: var = reduceLocalOffsets[node-numNodes];
665: } else {
666: var = reduceOffsets[node] - reduceFirstVar;
667: }
668: PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: u = %g v = %g", rank, node, PetscRealPart(reduceV[var]), PetscRealPart(reduceV[var+1]));
669: #ifdef DEBUG_ACC
670: PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: a = %g b = %g", rank, node, PetscRealPart(reduceA[var]), PetscRealPart(reduceA[var+1]));
671: #endif
672: } else {
673: if (node >= numNodes) {
674: var = localOffsets[node-numNodes];
675: } else {
676: var = offsets[node] - firstVar;
677: }
678: PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: u = %g v = %g", rank, node, PetscRealPart(v[var]), PetscRealPart(v[var+1]));
679: #ifdef DEBUG_ACC
680: PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: a = %g b = %g", rank, node, PetscRealPart(a[var]), PetscRealPart(a[var+1]));
681: #endif
682: }
683: MeshGetBoundaryNext(mesh, bd, PETSC_TRUE, &node);
684: }
685: }
686: VecRestoreArray(velGrid->ghostVec, &v);
687: #ifdef DEBUG_ACC
688: VecRestoreArray(accGrid->ghostVec, &a);
689: #endif
690: if (velGrid->reduceSystem == PETSC_TRUE) {
691: VecRestoreArray(velGrid->bdReduceVecCur, &reduceV);
692: VecRestoreArray(accGrid->bdReduceVecCur, &reduceA);
693: }
694: return(0);
695: }
697: static struct _MeshMoverOps MeshMoverOps = {/* Generic Operations */
698: MeshMoverSetup_Triangular_2D,
699: PETSC_NULL/* MeshMoverSetFromOptions */,
700: MeshMoverDestroy_Triangular_2D,
701: /* Mesh-Specific Operations */
702: MeshMoveMesh_Triangular_2D,
703: MeshUpdateNodeValues_Triangular_2D,
704: MeshReform_Triangular_2D};
706: EXTERN_C_BEGIN
707: int MeshMoverCreate_Triangular_2D(MeshMover mover, ParameterDict dict) {
711: PetscMemcpy(mover->ops, &MeshMoverOps, sizeof(struct _MeshMoverOps));
713: return(0);
714: }
715: EXTERN_C_END