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: #undef  __FUNCT__
 13: int MeshMoverSetup_Triangular_2D(MeshMover mover) {
 14:   Mesh               mesh       = mover->mesh;
 15:   Mesh               movingMesh = mesh;
 16:   int                dim        = mesh->dim;
 17:   int                numCorners = mesh->numCorners;
 18:   DiscretizationType dtype;
 19:   PetscTruth         opt;
 20:   int                ierr;

 23:   /* Determine order of interpolation */
 24:   if (numCorners == 3) {
 25:     dtype = DISCRETIZATION_TRIANGULAR_2D_LINEAR;
 26:   } else if (numCorners == 6) {
 27:     dtype = DISCRETIZATION_TRIANGULAR_2D_QUADRATIC;
 28:     ierr  = PetscOptionsHasName(mesh->prefix, "-mesh_move_linear", &opt);
 29:     if (opt == PETSC_TRUE) {
 30:       dtype = DISCRETIZATION_TRIANGULAR_2D_LINEAR;
 31:       /* Must play games here so that the coarse mesh does not get a MeshMover as well */
 32:       PetscObjectReference((PetscObject) mover);
 33:       MeshSetMover(mesh, PETSC_NULL);
 34:       MeshCoarsen(mesh, PETSC_NULL, &movingMesh);
 35:       MeshSetMover(mesh, mover);
 36:     }
 37:   } else {
 38:     SETERRQ1(PETSC_ERR_SUP, "Number of nodes per element %d not supported", numCorners);
 39:   }

 41:   /* Setup mesh velocity grid */
 42:   GridCreate(movingMesh, &mover->meshVelocityGrid);
 43:   GridAddField(mover->meshVelocityGrid, "Velocity", dtype, dim, PETSC_NULL);
 44:   GridAppendOptionsPrefix(mover->meshVelocityGrid, "mesh_vel_");
 45:   GridSetFromOptions(mover->meshVelocityGrid);
 46:   GridSetActiveField(mover->meshVelocityGrid, 0);
 47:   PetscLogObjectParent(movingMesh, mover->meshVelocityGrid);
 48:   /* ALE operators need boundary values */
 49:   mover->meshVelocityGrid->reduceElement = PETSC_TRUE;
 50:   switch(mover->meshVelocityMethod) {
 51:   case MESH_LAPLACIAN:
 52:     GridSetMatOperator(mover->meshVelocityGrid, 0, 0, LAPLACIAN,    1.0, PETSC_FALSE, PETSC_NULL);
 53:     break;
 54:   case MESH_WEIGHTED_LAPLACIAN:
 55:     GridSetMatOperator(mover->meshVelocityGrid, 0, 0, WEIGHTED_LAP, 1.0, PETSC_FALSE, PETSC_NULL);
 56:     break;
 57:   default:
 58:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid velocity calculation method %d", mover->meshVelocityMethod);
 59:   }

 61:   /* Setup mesh acceleration grid */
 62:   GridCreate(movingMesh, &mover->meshAccelerationGrid);
 63:   GridAddField(mover->meshAccelerationGrid, "Acceleration", dtype, dim, PETSC_NULL);
 64:   GridAppendOptionsPrefix(mover->meshAccelerationGrid, "mesh_accel_");
 65:   GridSetFromOptions(mover->meshAccelerationGrid);
 66:   GridSetActiveField(mover->meshAccelerationGrid, 0);
 67:   PetscLogObjectParent(mesh, mover->meshAccelerationGrid);
 68:   switch(mover->meshAccelerationMethod) {
 69:   case MESH_LAPLACIAN:
 70:     GridSetMatOperator(mover->meshAccelerationGrid, 0, 0, LAPLACIAN, 1.0, PETSC_FALSE, PETSC_NULL);
 71:     break;
 72:   case MESH_WEIGHTED_LAPLACIAN:
 73:     GridSetMatOperator(mover->meshVelocityGrid, 0, 0, WEIGHTED_LAP, 1.0, PETSC_FALSE, PETSC_NULL);
 74:     break;
 75:   default:
 76:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid acceleration calculation method%d", mover->meshAccelerationMethod);
 77:   }
 78:   /* This will cause the coarse mesh to be destroyed when no longer needed */
 79:   if (mesh != movingMesh) {
 80:     MeshDestroy(movingMesh);
 81:   }
 82:   return(0);
 83: }

 85: #undef  __FUNCT__
 87: int MeshMoverDestroy_Triangular_2D(MeshMover mover) {
 89:   return(0);
 90: }

 92: #undef  __FUNCT__
 94: int MeshMapVertices_Triangular_2D(GVec coarseVec, Mesh mesh, Vec vec) {
 95:   Partition    p   = mesh->part;
 96:   int          dim = mesh->dim;
 97:   Grid         coarseGrid;
 98:   Mesh         coarseMesh;
 99:   Partition    coarseP;
100:   AO           coarseMap;
101:   PetscTruth   reduceSystem;
102:   int          numCoarseNodes;
103:   int         *classes;
104:   int          firstVar;
105:   int         *offsets, *localOffsets;
106:   int          reduceFirstVar     = 0;
107:   int         *reduceOffsets      = PETSC_NULL;
108:   int         *reduceLocalOffsets = PETSC_NULL;
109:   int        **reduceFieldClasses = PETSC_NULL;
110:   PetscScalar *vals, *coarseVals, *reduceVals;
111:   int          node, cNode, nclass, coarseVar, j;
112:   int          ierr;

115:   GVecGetGrid(coarseVec,  &coarseGrid);
116:   GridGetMesh(coarseGrid, &coarseMesh);
117:   VecGetArray(vec,       &vals);
118:   VecGetArray(coarseVec, &coarseVals);
119:   coarseP              = coarseMesh->part;
120:   coarseMap            = coarseMesh->coarseMap;
121:   reduceSystem         = coarseGrid->reduceSystem;
122:   numCoarseNodes       = coarseMesh->numNodes;
123:   reduceFieldClasses   = coarseGrid->cm->reduceFieldClasses;
124:   classes              = coarseGrid->cm->classes;
125:   firstVar             = coarseGrid->order->firstVar[p->rank];
126:   offsets              = coarseGrid->order->offsets;
127:   localOffsets         = coarseGrid->order->localOffsets;
128:   if (reduceSystem == PETSC_TRUE) {
129:     reduceOffsets      = coarseGrid->reduceOrder->offsets;
130:     reduceLocalOffsets = coarseGrid->reduceOrder->localOffsets;
131:     reduceFirstVar     = coarseGrid->reduceOrder->firstVar[p->rank];
132:     VecGetArray(coarseGrid->bdReduceVecCur, &reduceVals);
133:   }
134:   for(cNode = 0; cNode < numCoarseNodes; cNode++) {
135:     /* Get corresponding node */
136:     PartitionLocalToGlobalNodeIndex(coarseP, cNode, &node);
137:     AOPetscToApplication(coarseMap, 1, &node);
138:     PartitionGlobalToLocalNodeIndex(p,       node,  &node);
139:     /* Calculate coarse offset */
140:     nclass = classes[cNode];
141:     if ((reduceSystem == PETSC_TRUE) && (reduceFieldClasses[0][nclass])) {
142:       if (cNode >= numCoarseNodes) {
143:         coarseVar = reduceLocalOffsets[cNode-numCoarseNodes];
144:       } else {
145:         coarseVar = reduceOffsets[cNode] - reduceFirstVar;
146:       }
147:       /* Set node value */
148:       for(j = 0; j < dim; j++) {
149:         vals[node*dim+j] = reduceVals[coarseVar+j];
150:       }
151:     } else {
152:       if (cNode >= numCoarseNodes) {
153:         coarseVar = localOffsets[cNode-numCoarseNodes];
154:       } else {
155:         coarseVar = offsets[cNode] - firstVar;
156:       }
157:       /* Set node value */
158:       for(j = 0; j < dim; j++) {
159:         vals[node*dim+j] = coarseVals[coarseVar+j];
160:       }
161:     }
162:   }
163:   VecRestoreArray(vec,       &vals);
164:   VecRestoreArray(coarseVec, &coarseVals);
165:   if (reduceSystem == PETSC_TRUE) {
166:     VecRestoreArray(coarseGrid->bdReduceVecCur, &reduceVals);
167:   }
168:   return(0);
169: }

171: #undef  __FUNCT__
173: int MeshInterpMidnodes_Triangular_2D(Mesh mesh, Vec vec) {
174:   Mesh_Triangular *tri         = (Mesh_Triangular *) mesh->data;
175:   int              dim         = mesh->dim;
176:   int              numCorners  = mesh->numCorners;
177:   int              numElements = mesh->numFaces;
178:   int             *elements    = tri->faces;
179:   PetscScalar     *vals;
180:   int              elem, corner, node, nNode1, nNode2, j;
181:   int              ierr;

184:   VecGetArray(vec, &vals);
185:   for(elem = 0; elem < numElements; elem++) {
186:     for(corner = 3; corner < numCorners; corner++) {
187:       node    = elements[elem*numCorners+corner];
188:       nNode1  = elements[elem*numCorners+((corner+1)%3)];
189:       nNode2  = elements[elem*numCorners+((corner+2)%3)];
190:       for(j = 0; j < dim; j++) {
191:         vals[node*dim+j] = 0.5*(vals[nNode1*dim+j] + vals[nNode2*dim+j]);
192:       }
193:     }
194:   }
195:   VecRestoreArray(vec, &vals);
196:   PetscLogFlops(2*dim*(numCorners-3)*numElements);
197:   return(0);
198: }

200: #undef  __FUNCT__
202: int MeshUpdateNodeValues_Triangular_2D(MeshMover mover, GVec sol, Vec vec, Vec ghostVec) {
203:   Mesh       mesh       = mover->mesh;
204:   Partition  p          = mesh->part;
205:   int        numCorners = mesh->numCorners;
206:   Grid       grid;
207:   int        numOverlapNodes, numFuncs;
208:   PetscTruth match;
209:   int        ierr;

212:   /* Check grid */
213:   GVecGetGrid(sol, &grid);
214:   PartitionGetNumOverlapNodes(p, &numOverlapNodes);
215:   DiscretizationGetNumFunctions(grid->fields[0].disc, &numFuncs);
216:   if (grid->reduceSystem == PETSC_TRUE) {
217:     if (numFuncs == numCorners) {
218:       if (grid->order->numOverlapVars + grid->reduceOrder->numOverlapVars != numOverlapNodes*2) {
219:         SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid number of mesh velocity unknowns %d should be %d",
220:                  grid->order->numOverlapVars + grid->reduceOrder->numOverlapVars, numOverlapNodes*2);
221:       }
222:     }
223:   } else {
224:     if (numFuncs == numCorners) {
225:       if (grid->order->numOverlapVars != numOverlapNodes*2) {
226:         SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid number of mesh velocity unknowns %d should be %d",
227:                  grid->order->numOverlapVars, numOverlapNodes*2);
228:       }
229:     }
230:   }
231:   /* Map coarse nodes to fine vertices */
232:   PetscTypeCompare((PetscObject) grid->fields[0].disc, DISCRETIZATION_TRIANGULAR_2D_LINEAR, &match);
233:   if ((numCorners == 6) && (match == PETSC_TRUE)) {
234:     MeshMapVertices_Triangular_2D(sol, mesh, vec);
235:   } else {
236:     VecCopy(sol, vec);
237:   }
238:   /* Copy to ghost vector */
239:   VecScatterBegin(vec, ghostVec, INSERT_VALUES, SCATTER_FORWARD, mover->meshVelocityScatter);
240:   VecScatterEnd(vec, ghostVec, INSERT_VALUES, SCATTER_FORWARD, mover->meshVelocityScatter);
241:   /* Interpolate velocities of midnodes */
242:   MeshInterpMidnodes_Triangular_2D(mesh, ghostVec);

244:   /* REMOVE THIS AFTER YOU REWRITE THE ALE STUFF */
245:   if (sol == mover->meshVelocitySol) {
246:     GridGlobalToLocal(mover->meshVelocityGrid, INSERT_VALUES, sol);
247:   } else if (sol == mover->meshAccelerationSol) {
248:     GridGlobalToLocal(mover->meshAccelerationGrid, INSERT_VALUES, sol);
249:   } else {
250:     SETERRQ(PETSC_ERR_SUP, "Unrecognized solution vector");
251:   }
252:   return(0);
253: }

255: #undef  __FUNCT__
257: int MeshUpdateCoords_Triangular_2D(Mesh mesh, PetscScalar *v, PetscScalar *a, double t) {
258:   Mesh_Triangular *tri   = (Mesh_Triangular *) mesh->data;
259:   Partition        p     = mesh->part;
260:   int              dim   = mesh->dim;
261:   double          *nodes = tri->nodes;
262:   int              numOverlapNodes;
263:   int              node;
264:   int              ierr;

267:   PartitionGetNumOverlapNodes(p, &numOverlapNodes);
268:   /* x' = x + v t + 0.5 a t^2 */
269:   if (mesh->isPeriodic == PETSC_TRUE) {
270:     for(node = 0; node < numOverlapNodes; node++) {
271:       nodes[node*2]   = MeshPeriodicX(mesh, nodes[node*2]   + PetscRealPart(v[node*2])*t   + 0.5*PetscRealPart(a[node*2])*t*t);
272:       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);
273: #ifdef DEBUG_BD_MOVEMENT
274:       if (tri->markers[node] > 0)
275:         PetscPrintf(PETSC_COMM_SELF, "node %d: u = %lf v = %lfn", node, v[node*2], v[node*2+1]);
276: #endif
277:     }
278:   } else {
279:     for(node = 0; node < numOverlapNodes*dim; node++) {
280:       nodes[node] += PetscRealPart(v[node])*t + 0.5*PetscRealPart(a[node])*t*t;
281: #ifdef DEBUG_BD_MOVEMENT
282:       if ((node%2 == 0) && (tri->markers[node/2] > 0))
283:         PetscPrintf(PETSC_COMM_SELF, "node %d: u = %lf v = %lfn", node/2, v[node], v[node+1]);
284: #endif
285:     }
286:   }
287:   PetscLogFlops(6*numOverlapNodes*dim);
288:   return(0);
289: }

291: #undef  __FUNCT__
293: int MeshSyncCoords_Triangular_2D(Mesh mesh, Mesh coarseMesh) {
294:   Mesh_Triangular *tri            = (Mesh_Triangular *) mesh->data;
295:   Mesh_Triangular *coarseTri      = (Mesh_Triangular *) coarseMesh->data;
296:   Partition        p              = mesh->part;
297:   Partition        coarseP        = coarseMesh->part;
298:   AO               coarseMap      = coarseMesh->coarseMap;
299:   int              dim            = mesh->dim;
300:   int              numCoarseNodes = coarseMesh->numNodes;
301:   double          *nodes          = tri->nodes;
302:   double          *coarseNodes    = coarseTri->nodes;
303:   int              node, cNode, j;
304:   int              ierr;

307:   for(cNode = 0; cNode < numCoarseNodes; cNode++) {
308:     PartitionLocalToGlobalNodeIndex(coarseP, cNode, &node);
309:     AOPetscToApplication(coarseMap, 1, &node);
310:     PartitionGlobalToLocalNodeIndex(p,       node,  &node);
311:     for(j = 0; j < dim; j++) {
312:       coarseNodes[cNode*dim+j] = nodes[node*dim+j];
313:     }
314:   }
315:   return(0);
316: }

318: #undef  __FUNCT__
320: int MeshMoveMesh_Triangular_2D(MeshMover mover, double t) {
321:   Mesh         mesh    = mover->mesh;
322:   Mesh         velMesh = mover->meshVelocityGrid->mesh;
323:   PetscScalar *v, *a;
324:   PetscDraw    draw;
325:   int          ierr;

328:   VecGetArray(mover->nodeVelocitiesGhost,    &v);
329:   VecGetArray(mover->nodeAccelerationsGhost, &a);

331:   /* Move mesh nodes */
332:   MeshUpdateCoords_Triangular_2D(mesh, v, a, t);
333:   /* Move mesh velocity and acceleration mesh nodes */
334:   MeshSyncCoords_Triangular_2D(mesh, velMesh);

336:   VecRestoreArray(mover->nodeVelocitiesGhost,    &v);
337:   VecRestoreArray(mover->nodeAccelerationsGhost, &a);
338: 
339:   /* Update bounding rectangle */
340:   MeshUpdateBoundingBox(mesh);

342:   /* View mesh */
343:   if (mover->movingMeshViewer) {
344:     PetscViewerDrawGetDraw(mover->movingMeshViewer, 0, &draw);
345:     PetscDrawClear(draw);
346:     PetscDrawSetTitle(draw, mover->movingMeshViewerCaption);
347:     MeshView(mesh, mover->movingMeshViewer);
348:     PetscViewerFlush(mover->movingMeshViewer);
349:     PetscDrawPause(draw);
350:   }
351:   return(0);
352: }

354: #undef  __FUNCT__
356: /*
357:   MeshReform_Triangular_2D - Reforms a two dimensional unstructured mesh
358:   using the original boundary.

360:   Input Parameters:
361: . mesh    - The mesh begin reformed
362: . refine  - Flag indicating whether to refine or recreate the mesh
363: . area    - A function which gives an area constraint when evaluated inside an element
364: . newBd   - Flag to determine whether the boundary should be generated (PETSC_TRUE) or taken from storage

366:   Output Parameter:
367: . newMesh - The reformed mesh

369: .keywords unstructured grid
370: .seealso GridReform(), MeshRefine(), MeshSetBoundary()
371: */
372: int MeshReform_Triangular_2D(Mesh mesh, PetscTruth refine, PointFunction area, PetscTruth newBd, Mesh *newMesh) {
373:   Mesh                     m;
374:   Mesh_Triangular         *tri = (Mesh_Triangular *) mesh->data;
375:   Partition                p   = mesh->part;
376:   Partition_Triangular_2D *q   = (Partition_Triangular_2D *) p->data;
377:   MeshBoundary2D           bdCtx;
378:   int                      numEdges     = mesh->numBdEdges;
379:   int                      origLocNodes = mesh->numBdNodes;
380:   int                      numLocNodes;
381:   int                     *firstBdEdge;
382:   int                     *firstBdNode;
383:   int                     *firstHole;
384:   int                     *numRecvNodes;
385:   int                     *numRecvMarkers;
386:   int                     *numRecvSegments;
387:   int                     *numRecvSegMarkers;
388:   int                     *numRecvHoles;
389:   int                     *nodeOffsets;
390:   int                     *markerOffsets;
391:   int                     *segmentOffsets;
392:   int                     *segMarkerOffsets;
393:   int                     *holeOffsets;
394:   double                  *locNodes;
395:   int                     *locMarkers;
396:   int                     *locSegments;
397:   int                     *locSegMarkers;
398:   double                  *nodes;
399:   int                     *markers;
400:   int                     *segments;
401:   int                     *segMarkers;
402:   double                  *holes;
403:   int                     *globalNodes;
404:   int                     *aOrder;
405:   int                      rank, numProcs;
406:   int                      numLocSegments, numSegments, numNodes, numHoles;
407:   int                      proc, bd, gNode, lNode, bdNode, edge, bdEdge, edgeEnd;
408:   int                      ierr;

411:   MPI_Comm_size(mesh->comm, &numProcs);
412:   MPI_Comm_rank(mesh->comm, &rank);
413:   if ((newBd == PETSC_TRUE) || (mesh->bdCtxNew == PETSC_NULL)) {
414:     MPI_Reduce(&mesh->numHoles, &numHoles, 1, MPI_INT, MPI_SUM, 0, mesh->comm);
415:     PetscMalloc((numProcs+1)   * sizeof(int),    &firstBdEdge);
416:     PetscMalloc((numProcs+1)   * sizeof(int),    &firstBdNode);
417:     PetscMalloc((numProcs+1)   * sizeof(int),    &firstHole);
418:     PetscMalloc(numProcs       * sizeof(int),    &numRecvNodes);
419:     PetscMalloc(numProcs       * sizeof(int),    &numRecvMarkers);
420:     PetscMalloc(numProcs       * sizeof(int),    &numRecvSegments);
421:     PetscMalloc(numProcs       * sizeof(int),    &numRecvSegMarkers);
422:     PetscMalloc(numProcs       * sizeof(int),    &numRecvHoles);
423:     PetscMalloc(numProcs       * sizeof(int),    &nodeOffsets);
424:     PetscMalloc(numProcs       * sizeof(int),    &markerOffsets);
425:     PetscMalloc(numProcs       * sizeof(int),    &segmentOffsets);
426:     PetscMalloc(numProcs       * sizeof(int),    &segMarkerOffsets);
427:     PetscMalloc(numProcs       * sizeof(int),    &holeOffsets);
428:     PetscMalloc(origLocNodes*2 * sizeof(double), &locNodes);
429:     PetscMalloc(origLocNodes   * sizeof(int),    &locMarkers);
430:     PetscMalloc(numEdges*2     * sizeof(int),    &locSegments);
431:     PetscMalloc(numEdges       * sizeof(int),    &locSegMarkers);
432:     PetscMalloc(origLocNodes   * sizeof(int),    &globalNodes);

434:     /* Make all node numbers global */
435:     numLocNodes    = 0;
436:     numLocSegments = 0;
437:     for(bd = 0; bd < mesh->numBd; bd++) {
438:       /* Get boundary nodes */
439:       for(bdNode = tri->bdBegin[bd]; bdNode < tri->bdBegin[bd+1]; bdNode++) {
440:         gNode = tri->bdNodes[bdNode];
441:         lNode = gNode - q->firstNode[p->rank];
442:         if ((lNode >= 0) && (lNode < mesh->numNodes)) {
443:           if (tri->degrees[lNode] > 0) {
444:             locNodes[numLocNodes*2]   = tri->nodes[lNode*2];
445:             locNodes[numLocNodes*2+1] = tri->nodes[lNode*2+1];
446:             locMarkers[numLocNodes]   = tri->bdMarkers[bd];
447:             globalNodes[numLocNodes]  = gNode;
448:             numLocNodes++;
449:           }
450:         }
451:       }

453:       /* Get boundary segments */
454:       for(bdEdge = tri->bdEdgeBegin[bd]; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++) {
455:         edge = tri->bdEdges[bdEdge] - q->firstEdge[rank];
456:         if ((edge >= 0) && (edge < mesh->numEdges)) {
457:           for(edgeEnd = 0; edgeEnd < 2; edgeEnd++) {
458:             lNode = tri->edges[edge*2+edgeEnd];
459:             ierr  = PartitionLocalToGlobalNodeIndex(p, lNode, &gNode);
460:             locSegments[numLocSegments*2+edgeEnd] = gNode;
461:           }
462:           locSegMarkers[numLocSegments]   = tri->bdMarkers[bd];
463:           numLocSegments++;
464:         }
465:       }
466:     }

468:     /* Assemble offsets */
469:     MPI_Allgather(&numLocSegments, 1, MPI_INT, &firstBdEdge[1], 1, MPI_INT, p->comm);
470:     firstBdEdge[0] = 0;
471:     for(proc = 1; proc <= numProcs; proc++)
472:       firstBdEdge[proc] = firstBdEdge[proc] + firstBdEdge[proc-1];
473:     numSegments = firstBdEdge[numProcs];
474:     if ((rank == 0) && (numSegments != numEdges)) SETERRQ(PETSC_ERR_PLIB, "Invalid number of boundary segments");
475:     MPI_Allgather(&numLocNodes,    1, MPI_INT, &firstBdNode[1], 1, MPI_INT, p->comm);
476:     firstBdNode[0] = 0;
477:     for(proc = 1; proc <= numProcs; proc++)
478:       firstBdNode[proc] = firstBdNode[proc] + firstBdNode[proc-1];
479:     numNodes = firstBdNode[numProcs];
480:     MPI_Allgather(&mesh->numHoles,  1, MPI_INT, &firstHole[1],   1, MPI_INT, p->comm);
481:     firstHole[0] = 0;
482:     for(proc = 1; proc <= numProcs; proc++)
483:       firstHole[proc] = firstHole[proc] + firstHole[proc-1];
484:     if ((rank == 0) && (numHoles != firstHole[numProcs])) SETERRQ(PETSC_ERR_PLIB, "Invalid number of holes");

486:     for(proc = 0; proc < numProcs; proc++) {
487:       numRecvNodes[proc]      = (firstBdNode[proc+1] - firstBdNode[proc])*2;
488:       numRecvMarkers[proc]    = (firstBdNode[proc+1] - firstBdNode[proc]);
489:       numRecvSegments[proc]   = (firstBdEdge[proc+1] - firstBdEdge[proc])*2;
490:       numRecvSegMarkers[proc] = (firstBdEdge[proc+1] - firstBdEdge[proc]);
491:       numRecvHoles[proc]      = (firstHole[proc+1]   - firstHole[proc])*2;
492:     }
493:     nodeOffsets[0]      = 0;
494:     markerOffsets[0]    = 0;
495:     segmentOffsets[0]   = 0;
496:     segMarkerOffsets[0] = 0;
497:     holeOffsets[0]      = 0;
498:     for(proc = 1; proc < numProcs; proc++) {
499:       nodeOffsets[proc]      = numRecvNodes[proc-1]      + nodeOffsets[proc-1];
500:       markerOffsets[proc]    = numRecvMarkers[proc-1]    + markerOffsets[proc-1];
501:       segmentOffsets[proc]   = numRecvSegments[proc-1]   + segmentOffsets[proc-1];
502:       segMarkerOffsets[proc] = numRecvSegMarkers[proc-1] + segMarkerOffsets[proc-1];
503:       holeOffsets[proc]      = numRecvHoles[proc-1]      + holeOffsets[proc-1];
504:     }

506:     if (rank == 0) {
507:       PetscMalloc(numNodes*2 * sizeof(double), &nodes);
508:       PetscMalloc(numNodes   * sizeof(int),    &markers);
509:       PetscMalloc(numEdges*2 * sizeof(int),    &segments);
510:       PetscMalloc(numEdges   * sizeof(int),    &segMarkers);
511:       if (numHoles) {
512:         ierr  = PetscMalloc(numHoles*2 * sizeof(double), &holes);
513:       } else {
514:         holes = PETSC_NULL;
515:       }
516:     }

518:     /* Make new numbering for boundary nodes */
519:     PetscMalloc(numNodes * sizeof(int), &aOrder);
520:     MPI_Allgatherv(globalNodes, numLocNodes, MPI_INT, aOrder, numRecvMarkers, markerOffsets, MPI_INT, p->comm);
521: 

523:     /* Renumber segments */
524:     for(edge = 0; edge < numLocSegments*2; edge++) {
525:       gNode = locSegments[edge];
526:       for(bdNode = 0; bdNode < numNodes; bdNode++)
527:         if (aOrder[bdNode] == gNode)
528:           break;
529:       if (bdNode == numNodes) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary node");
530:       locSegments[edge] = bdNode;
531:     }
532:     PetscFree(aOrder);

534:     MPI_Gatherv(locNodes,      numLocNodes*2,    MPI_DOUBLE, nodes,      numRecvNodes,      nodeOffsets,      MPI_DOUBLE, 0, mesh->comm);
535: 
536:     MPI_Gatherv(locMarkers,    numLocNodes,      MPI_INT,    markers,    numRecvMarkers,    markerOffsets,    MPI_INT,    0, mesh->comm);
537: 
538:     MPI_Gatherv(locSegments,   numLocSegments*2, MPI_INT,    segments,   numRecvSegments,   segmentOffsets,   MPI_INT,    0, mesh->comm);
539: 
540:     MPI_Gatherv(locSegMarkers, numLocSegments,   MPI_INT,    segMarkers, numRecvSegMarkers, segMarkerOffsets, MPI_INT,    0, mesh->comm);
541: 
542:     MPI_Gatherv(mesh->holes,   mesh->numHoles*2, MPI_DOUBLE, holes,      numRecvHoles,      holeOffsets,      MPI_DOUBLE, 0, mesh->comm);
543: 

545:     /* Cleanup */
546:     PetscFree(firstBdEdge);
547:     PetscFree(firstBdNode);
548:     PetscFree(firstHole);
549:     PetscFree(numRecvNodes);
550:     PetscFree(numRecvMarkers);
551:     PetscFree(numRecvSegments);
552:     PetscFree(numRecvSegMarkers);
553:     PetscFree(numRecvHoles);
554:     PetscFree(nodeOffsets);
555:     PetscFree(markerOffsets);
556:     PetscFree(segmentOffsets);
557:     PetscFree(segMarkerOffsets);
558:     PetscFree(holeOffsets);
559:     PetscFree(locNodes);
560:     PetscFree(locMarkers);
561:     PetscFree(locSegments);
562:     PetscFree(locSegMarkers);
563:     PetscFree(globalNodes);

565:     bdCtx.numVertices = numNodes;
566:     bdCtx.vertices    = nodes;
567:     bdCtx.markers     = markers;
568:     bdCtx.numSegments = numEdges;
569:     bdCtx.segments    = segments;
570:     bdCtx.segMarkers  = segMarkers;
571:     bdCtx.numHoles    = numHoles;
572:     bdCtx.holes       = holes;
573:   } else {
574:     bdCtx.numVertices = mesh->bdCtxNew->numVertices;
575:     bdCtx.vertices    = mesh->bdCtxNew->vertices;
576:     bdCtx.markers     = mesh->bdCtxNew->markers;
577:     bdCtx.numSegments = mesh->bdCtxNew->numSegments;
578:     bdCtx.segments    = mesh->bdCtxNew->segments;
579:     bdCtx.segMarkers  = mesh->bdCtxNew->segMarkers;
580:     bdCtx.numHoles    = mesh->bdCtxNew->numHoles;
581:     bdCtx.holes       = mesh->bdCtxNew->holes;
582:   }
583:   bdCtx.numBd = mesh->numBd;

585:   MeshCreate(mesh->comm, &m);
586:   MeshSetDimension(m, 2);
587:   MeshSetPeriodicDimension(m, 0, mesh->isPeriodicDim[0]);
588:   MeshSetPeriodicDimension(m, 1, mesh->isPeriodicDim[1]);
589:   MeshSetNumCorners(m, mesh->numCorners);
590:   MeshSetBoundary(m, &bdCtx);
591:   MeshSetType(m, mesh->type_name);
592:   PetscObjectSetName((PetscObject) m, "Mesh");
593:   *newMesh = m;

595:   if (rank == 0) {
596:     if ((newBd == PETSC_TRUE) || (mesh->bdCtxNew == PETSC_NULL)) {
597:       PetscFree(nodes);
598:       PetscFree(markers);
599:       PetscFree(segments);
600:       PetscFree(segMarkers);
601:       if (numHoles) {
602:         PetscFree(holes);
603:       }
604:     } else {
605: #if 0
606:       MeshBoundaryDestroy(mesh->bdCtxNew);
607: #endif
608:       mesh->bdCtxNew = PETSC_NULL;
609:     }
610:   }

612: #ifdef PETSC_USE_BOPT_g
613:   PetscTrValid(__LINE__, __FUNCT__, __FILE__, __SDIR__);
614: #endif
615:   return(0);
616: }

618: #undef  __FUNCT__
620: int MeshDebugBoundary_Triangular_2D_Private(Mesh mesh, GVec nodeVelocities) {
621:   Mesh_Triangular *tri      = (Mesh_Triangular *) mesh->data;
622:   int              numNodes = mesh->numNodes;
623:   double          *nodes    = tri->nodes;
624:   int             *markers  = tri->markers;
625:   int              rank     = mesh->part->rank;
626:   PetscScalar     *vel;
627:   int              node;
628:   int              ierr;

631:   VecGetArray(nodeVelocities, &vel);
632:   for(node = 0; node < numNodes; node++) {
633:     if (markers[node] > 0) {
634:       PetscPrintf(PETSC_COMM_SELF, "[%d]node %d at (%g,%g) vel: (%g,%g)n",
635:                   rank, node, nodes[node*2], nodes[node*2+1], PetscRealPart(vel[node*2]), PetscRealPart(vel[node*2+1]));
636:     }
637:   }
638:   VecRestoreArray(nodeVelocities, &vel);
639:   return(0);
640: }

642: #undef  __FUNCT__
644: int MeshDebugBoundary_Triangular_2D(MeshMover mover) {
645:   Mesh             mesh               = mover->mesh;
646:   Mesh_Triangular *tri                = (Mesh_Triangular *) mesh->data;
647:   Grid             velGrid            = mover->meshVelocityGrid;
648:   Grid             accGrid            = mover->meshAccelerationGrid;
649:   int              numBd              = mesh->numBd;
650:   int             *bdMarkers          = tri->bdMarkers;
651:   int              numNodes           = mesh->numNodes;
652:   int              rank               = mesh->part->rank;
653:   int            **reduceFieldClasses = velGrid->cm->reduceFieldClasses;
654:   int             *classes            = velGrid->cm->classes;
655:   int              firstVar           = velGrid->order->firstVar[rank];
656:   int             *offsets            = velGrid->order->offsets;
657:   int             *localOffsets       = velGrid->order->localOffsets;
658:   int              reduceFirstVar     = 0;
659:   int             *reduceOffsets      = PETSC_NULL;
660:   int             *reduceLocalOffsets = PETSC_NULL;
661:   PetscScalar     *v, *reduceV, *reduceA;
662:   int              b, bd, node, nclass, var;
663:   int              ierr;

666:   VecGetArray(velGrid->ghostVec, &v);
667: #ifdef DEBUG_ACC
668:   VecGetArray(accGrid->ghostVec, &a);
669: #endif
670:   if (velGrid->reduceSystem == PETSC_TRUE) {
671:     reduceOffsets      = velGrid->reduceOrder->offsets;
672:     reduceLocalOffsets = velGrid->reduceOrder->localOffsets;
673:     reduceFirstVar     = velGrid->reduceOrder->firstVar[rank];
674:     VecGetArray(velGrid->bdReduceVecCur, &reduceV);
675:     VecGetArray(accGrid->bdReduceVecCur, &reduceA);
676:   }
677:   for(b = 0; b < numBd; b++) {
678:     bd = bdMarkers[b];
679:     if (bd < 0) continue;

681:     MeshGetBoundaryStart(mesh, bd, PETSC_TRUE, &node);
682:     while (node >= 0) {
683:       /* Print the value of velocity on boundary node */
684:       nclass = classes[node];
685:       if ((reduceFieldClasses != PETSC_NULL) && (reduceFieldClasses[0][nclass])) {
686:         if (node >= numNodes) {
687:           var = reduceLocalOffsets[node-numNodes];
688:         } else {
689:           var = reduceOffsets[node] - reduceFirstVar;
690:         }
691:         PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: u = %g v = %g", rank, node, PetscRealPart(reduceV[var]), PetscRealPart(reduceV[var+1]));
692: #ifdef DEBUG_ACC
693:         PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: a = %g b = %g", rank, node, PetscRealPart(reduceA[var]), PetscRealPart(reduceA[var+1]));
694: #endif
695:       } else {
696:         if (node >= numNodes) {
697:           var = localOffsets[node-numNodes];
698:         } else {
699:           var = offsets[node] - firstVar;
700:         }
701:         PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: u = %g v = %g", rank, node, PetscRealPart(v[var]), PetscRealPart(v[var+1]));
702: #ifdef DEBUG_ACC
703:         PetscPrintf(PETSC_COMM_SELF, "[%d]node %d: a = %g b = %g", rank, node, PetscRealPart(a[var]), PetscRealPart(a[var+1]));
704: #endif
705:       }
706:       MeshGetBoundaryNext(mesh, bd, PETSC_TRUE, &node);
707:     }
708:   }
709:   VecRestoreArray(velGrid->ghostVec, &v);
710: #ifdef DEBUG_ACC
711:   VecRestoreArray(accGrid->ghostVec, &a);
712: #endif
713:   if (velGrid->reduceSystem == PETSC_TRUE) {
714:     VecRestoreArray(velGrid->bdReduceVecCur, &reduceV);
715:     VecRestoreArray(accGrid->bdReduceVecCur, &reduceA);
716:   }
717:   return(0);
718: }

720: static struct _MeshMoverOps MeshMoverOps = {/* Generic Operations */
721:                                             MeshMoverSetup_Triangular_2D,
722:                                             PETSC_NULL/* MeshMoverSetFromOptions */,
723:                                             MeshMoverDestroy_Triangular_2D,
724:                                             /* Mesh-Specific Operations */
725:                                             MeshMoveMesh_Triangular_2D,
726:                                             MeshUpdateNodeValues_Triangular_2D,
727:                                             MeshReform_Triangular_2D};

729: EXTERN_C_BEGIN
730: #undef  __FUNCT__
732: int MeshMoverCreate_Triangular_2D(MeshMover mover, ParameterDict dict) {

736:   PetscMemcpy(mover->ops, &MeshMoverOps, sizeof(struct _MeshMoverOps));

738:   return(0);
739: }
740: EXTERN_C_END