Actual source code: tri2d_Triangle.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: tri2d_Triangle.c,v 1.6 2000/10/17 13:48:57 knepley Exp $";
  3: #endif

  5: #include "src/mesh/impls/triangular/2d/2dimpl.h"         /*I "mesh.h" I*/
  6: #include "tri2d_Triangle.h"

  8: #ifdef PETSC_HAVE_TRIANGLE

 10: #define ANSI_DECLARATORS
 11: #include "triangle.h"

 13: extern int MeshSetupBoundary_Triangular_2D(Mesh);
 14: extern int MeshDebug_Triangular_2D(Mesh, PetscTruth);
 15: extern int MeshSetupSupport_Triangular_2D(Mesh);
 16: extern int MeshAssemble_Private(Mesh, double **, int **, int **, int **);

 18: /*-------------------------------------------- Standard Interface ---------------------------------------------------*/
 19: /* MeshInitInput_Triangle
 20:   This function initializes the input structure for Triangle.

 22:   Collective on Mesh

 24:   Output Parameter:
 25: . inputCtx - Structure for communicating with Triangle

 27:   Level: developer

 29: .keywords mesh, Triangle
 30: .seealso MeshInitOutput_Triangle()
 31: */
 32: int MeshInitInput_Triangle(struct triangulateio *inputCtx)
 33: {
 35:   /* Setup the input structure */
 36:   inputCtx->numberofpoints          = 0;
 37:   inputCtx->numberofpointattributes = 0;
 38:   inputCtx->pointlist               = PETSC_NULL;
 39:   inputCtx->pointattributelist      = PETSC_NULL;
 40:   inputCtx->pointmarkerlist         = PETSC_NULL;
 41:   inputCtx->numberofsegments        = 0;
 42:   inputCtx->segmentlist             = PETSC_NULL;
 43:   inputCtx->segmentmarkerlist       = PETSC_NULL;
 44:   inputCtx->numberofholes           = 0;
 45:   inputCtx->holelist                = PETSC_NULL;
 46:   inputCtx->numberofregions         = 0;
 47:   inputCtx->regionlist              = PETSC_NULL;
 48:   return(0);
 49: }

 51: /* MeshInitRefineInput_Triangle
 52:   This function fills the input structure for Triangle with mesh information before refinement.

 54:   Collective on Mesh

 56:   Input Parameters:
 57: + mesh    - The mesh being refined
 58: - area     - A function which gives an area constraint when evaluated inside an element

 60:   Output Parameter:
 61: . inputCtx - Structure for communicating with Triangle

 63:   Level: developer

 65: .keywords mesh, Triangle
 66: .seealso MeshInitInput_Triangle(), MeshInitOutput_Triangle()
 67: */
 68: int MeshInitRefineInput_Triangle(Mesh mesh, PointFunction area, struct triangulateio *inputCtx)
 69: {
 70:   Mesh_Triangular         *tri        = (Mesh_Triangular *) mesh->data;
 71:   Partition                p          = mesh->part;
 72:   Partition_Triangular_2D *q          = (Partition_Triangular_2D *) p->data;
 73:   int                      numNodes   = q->numNodes;
 74:   int                      numFaces   = p->numElements;
 75:   int                      numCorners = mesh->numCorners;
 76:   int                      rank       = p->rank;
 77:   double                  *nodes;
 78:   int                     *markers, *faces, *edges;
 79:   int                     *segments, *segmarkers;
 80:   double                   x, y;                    /* Centroid of the current element */
 81:   int                      elem, bd, edge;
 82:   int                      ierr;

 85:   MeshAssemble_Private(mesh, &nodes, &markers, &faces, &edges);

 87:   if (rank == 0) {
 88:     PetscMalloc(mesh->numBdEdges*2 * sizeof(int), &segments);
 89:     PetscMalloc(mesh->numBdEdges   * sizeof(int), &segmarkers);

 91:     /* Needed to maintain boundary markers on edges */
 92:     for(bd = 0; bd < mesh->numBd; bd++) {
 93:       for(edge = tri->bdEdgeBegin[bd]; edge < tri->bdEdgeBegin[bd+1]; edge++) {
 94:         /* Make all node numbers global */
 95:         segments[edge*2]   = edges[tri->bdEdges[edge]*2];
 96:         segments[edge*2+1] = edges[tri->bdEdges[edge]*2+1];
 97:         segmarkers[edge]   = tri->bdMarkers[bd];
 98:       }
 99:     }
100:     PetscFree(edges);

102:     if (q->nodeOrdering != PETSC_NULL) {
103:       /* We must globally renumber and permute so that midnodes come after vertices */
104:       AOPetscToApplication(q->nodeOrdering, numFaces*numCorners, faces);
105:       AOPetscToApplicationPermuteReal(q->nodeOrdering, 2, nodes);
106:       AOPetscToApplicationPermuteInt(q->nodeOrdering,    1, markers);
107:       /* We can globally renumber only the segments, not the whole edge array */
108:       AOPetscToApplication(q->nodeOrdering, mesh->numBdEdges*2, segments);
109:     }
110:     if (mesh->nodeOrdering != PETSC_NULL) {
111:       AOPetscToApplication(mesh->nodeOrdering, numFaces*numCorners, faces);
112:       AOPetscToApplicationPermuteReal(mesh->nodeOrdering, 2, nodes);
113:       AOPetscToApplicationPermuteInt(mesh->nodeOrdering,    1, markers);
114:       AOPetscToApplication(mesh->nodeOrdering, mesh->numBdEdges*2, segments);
115:     }

117:     /* Setup the input structure */
118:     inputCtx->numberofpoints             = numNodes;
119:     inputCtx->numberofpointattributes    = 0;
120:     inputCtx->pointlist                  = nodes;
121:     inputCtx->pointattributelist         = PETSC_NULL;
122:     inputCtx->pointmarkerlist            = markers;
123:     inputCtx->numberoftriangles          = numFaces;
124:     inputCtx->numberofcorners            = numCorners;
125:     inputCtx->trianglelist               = faces;
126:     inputCtx->numberoftriangleattributes = 0;
127:     inputCtx->triangleattributelist      = PETSC_NULL;
128:     inputCtx->trianglearealist           = PETSC_NULL;
129:     inputCtx->numberofsegments           = mesh->numBdEdges;
130:     inputCtx->segmentlist                = segments;
131:     inputCtx->segmentmarkerlist          = segmarkers;
132:     inputCtx->numberofholes              = 0;
133:     inputCtx->holelist                   = PETSC_NULL;
134:     inputCtx->numberofregions            = 0;
135:     inputCtx->regionlist                 = PETSC_NULL;

137:     /* Create list of area constraints */
138:     PetscMalloc(numFaces * sizeof(double), &inputCtx->trianglearealist);
139:     for(elem = 0; elem < numFaces; elem++) {
140:       PetscReal maxArea;

142:       /* Calculate the centroid of the element */
143:       x = MeshPeriodicX(mesh, (nodes[faces[elem*numCorners]*2]   + nodes[faces[elem*numCorners+1]*2]   +
144:                         nodes[faces[elem*numCorners+2]*2])  /3.0);
145:       y = MeshPeriodicY(mesh, (nodes[faces[elem*numCorners]*2+1] + nodes[faces[elem*numCorners+1]*2+1] +
146:                         nodes[faces[elem*numCorners+2]*2+1])/3.0);
147:       (*area)(1, 1, &x, &y, PETSC_NULL, &maxArea, mesh->usr);
148:       inputCtx->trianglearealist[elem] = maxArea;
149:     }
150:   }
151: #ifdef PETSC_USE_BOPT_g
152: #endif

154:   return(0);
155: }

157: /* MeshFinalizeRefineInput_Triangle
158:   This function cleans up the input structure for Triangle after refinement.

160:   Collective on Mesh

162:   Input Parameters:
163: . mesh     - The mesh being refined

165:   Output Parameters:
166: . inputCtx - Structure for communicating with Triangle

168:   Level: developer

170: .keywords mesh, Triangle
171: .seealso MeshInitRefineInput_Triangle(), MeshInitInput_Triangle()
172: */
173: int MeshFinalizeRefineInput_Triangle(Mesh mesh, struct triangulateio *inputCtx)
174: {

178:   if (mesh->part->rank == 0) {
179:     PetscFree(inputCtx->pointlist);
180:     PetscFree(inputCtx->pointmarkerlist);
181:     PetscFree(inputCtx->trianglelist);
182:     PetscFree(inputCtx->segmentlist);
183:     PetscFree(inputCtx->segmentmarkerlist);
184:     PetscFree(inputCtx->trianglearealist);
185:   }
186:   return(0);
187: }

189: /* MeshInitOutput_Triangle
190:   This function initializes the output structure for Triangle.

192:   Collective on Mesh

194:   Output Parameter:
195: . outputCtx - Structure for communicating with Triangle

197:   Level: developer

199: .keywords mesh, Triangle
200: .seealso MeshInitInput_Triangle()
201: */
202: int MeshInitOutput_Triangle(struct triangulateio *outputCtx)
203: {
205:   /* Setup the output structure */
206:   outputCtx->pointlist             = PETSC_NULL;
207:   outputCtx->pointattributelist    = PETSC_NULL;
208:   outputCtx->pointmarkerlist       = PETSC_NULL;
209:   outputCtx->trianglelist          = PETSC_NULL;
210:   outputCtx->triangleattributelist = PETSC_NULL;
211:   outputCtx->neighborlist          = PETSC_NULL;
212:   outputCtx->segmentlist           = PETSC_NULL;
213:   outputCtx->segmentmarkerlist     = PETSC_NULL;
214:   outputCtx->edgelist              = PETSC_NULL;
215:   outputCtx->edgemarkerlist        = PETSC_NULL;
216:   return(0);
217: }

219: /* MeshDistribute_Triangle
220:   This function distributes the output from Triangle identically among processors.

222:   Collective on Mesh

224:   Output Parameters:
225: + mesh - The mesh
226: - out  - Structure for communicating with Triangle

228:   Level: developer

230: .keywords mesh, Triangle
231: .seealso MeshInitInput_Triangle
232: */
233: int MeshDistribute_Triangle(Mesh mesh, struct triangulateio *out)
234: {
235:   int rank;

239:   /* Communicate values */
240:   MPI_Comm_rank(mesh->comm, &rank);
241:   MPI_Bcast(&out->numberofpoints,    1, MPI_INT, 0, mesh->comm);
242:   MPI_Bcast(&out->numberofedges,     1, MPI_INT, 0, mesh->comm);
243:   MPI_Bcast(&out->numberoftriangles, 1, MPI_INT, 0, mesh->comm);
244:   MPI_Bcast(&out->numberofcorners,   1, MPI_INT, 0, mesh->comm);
245:   if (rank != 0) {
246:     out->numberofholes   = 0;
247:     out->holelist        = PETSC_NULL;
248:     PetscMalloc(out->numberofpoints*2                       * sizeof(double), &out->pointlist);
249:     PetscMalloc(out->numberofpoints                         * sizeof(int),    &out->pointmarkerlist);
250:     PetscMalloc(out->numberofedges*2                        * sizeof(int),    &out->edgelist);
251:     PetscMalloc(out->numberofedges                          * sizeof(int),    &out->edgemarkerlist);
252:     PetscMalloc(out->numberoftriangles*3                    * sizeof(int),    &out->neighborlist);
253:     PetscMalloc(out->numberoftriangles*out->numberofcorners * sizeof(int),    &out->trianglelist);
254:   }
255:   MPI_Bcast(out->pointlist,       out->numberofpoints*2,    MPI_DOUBLE, 0, mesh->comm);
256:   MPI_Bcast(out->pointmarkerlist, out->numberofpoints,      MPI_INT,    0, mesh->comm);
257:   MPI_Bcast(out->edgelist,        out->numberofedges*2,     MPI_INT,    0, mesh->comm);
258:   MPI_Bcast(out->edgemarkerlist,  out->numberofedges,       MPI_INT,    0, mesh->comm);
259:   MPI_Bcast(out->neighborlist,    out->numberoftriangles*3, MPI_INT,    0, mesh->comm);
260:   MPI_Bcast(out->trianglelist,    out->numberoftriangles*out->numberofcorners, MPI_INT, 0, mesh->comm);
261:   return(0);
262: }

264: /* MeshCreate_Triangle
265:   This function creates a 2D unstructured mesh using Triangle.

267:   Collective on Mesh

269:   Input Parameter:
270: . numCorners  - The number of nodes in an element

272:   Input Parameters from bdCtx:
273: + numBD       - The number of closed boundaries in the geometry, or different markers
274: . numVertices - The number of boundary points
275: . vertices    - The (x,y) coordinates of the boundary points
276: . markers     - The boundary markers for nodes, 0 indicates an interior point, each boundary must have a different marker
277: . numSegments - The number of boundary segments
278: . segments    - The endpoints of boundary segments or PETSC_NULL
279: . segMarkers  - The boundary markers for each segment
280: . numHoles    - The number of holes
281: - holes       - The (x,y) coordinates of holes or PETSC_NULL

283:   Output Parameter:
284: . mesh        - The new mesh created by Triangle

286:   Level: developer

288: .keywords mesh, Triangle
289: .seealso MeshInitInput_Triangle
290: */
291: int MeshCreate_Triangle(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
292: {
293:   Mesh_Triangular     *tri = (Mesh_Triangular *) mesh->data;
294:   struct triangulateio in;            /* Input  for Triangle mesh generator */
295:   struct triangulateio out;           /* Output for Triangle mesh generator */
296:   char                 cmd_line[256]; /* Command line for Triangle */
297:   int                  rank;
298:   PetscTruth           opt;
299:   int                  ierr;

302:   /* Initialize communication structures for Triangle */
303:   MeshInitInput_Triangle(&in);
304:   MeshInitOutput_Triangle(&out);

306:   MPI_Comm_rank(mesh->comm, &rank);
307:   if (rank == 0) {
308:     /* Setup boundaries and holes */
309:     in.numberofpoints      = bdCtx->numVertices;
310:     if (bdCtx->numVertices > 0) {
313:       in.pointlist         = bdCtx->vertices;
314:       in.pointmarkerlist   = bdCtx->markers;
315:     }
316:     in.numberofsegments    = bdCtx->numSegments;
317:     if (bdCtx->numSegments > 0) {
320:       in.segmentlist       = bdCtx->segments;
321:       in.segmentmarkerlist = bdCtx->segMarkers;
322:     }
323:     in.numberofholes       = bdCtx->numHoles;
324:     if (bdCtx->numHoles    > 0) {
325:       /* We keep these */
327:       PetscMalloc(in.numberofholes*2 * sizeof(double), &in.holelist);
328:       PetscMemcpy(in.holelist, bdCtx->holes, in.numberofholes*2 * sizeof(double));
329:     }

331:     /* Generate the inital coarse mesh and check whether we should create midnodes */
332:     if (numCorners == 3) {
333:       PetscStrcpy(cmd_line, "pqcenzQ");
334:     } else if (numCorners == 6) {
335:       PetscStrcpy(cmd_line, "pqcenzo2Q");
336:     } else {
337:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", numCorners);
338:     }
339:     PetscOptionsGetString(PETSC_NULL, "-triangle_cmd", cmd_line, 255, &opt);
340:     triangulate(cmd_line, &in, &out, PETSC_NULL);

342:     /* Get rid of extra information */
343:     PetscFree(out.segmentlist);
344:     PetscFree(out.segmentmarkerlist);
345:   }

347:   /* Communicate values */
348:   MeshDistribute_Triangle(mesh, &out);

350:   /* Store the information */
351:   mesh->numBd       = bdCtx->numBd;
352:   mesh->numNodes    = out.numberofpoints;
353:   mesh->numEdges    = out.numberofedges;
354:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
355:   mesh->numFaces    = out.numberoftriangles;
356:   mesh->numCorners  = out.numberofcorners;
357:   mesh->numHoles    = out.numberofholes;
358:   tri->nodes       = out.pointlist;
359:   tri->markers     = out.pointmarkerlist;
360:   tri->edges       = out.edgelist;
361:   tri->edgemarkers = out.edgemarkerlist;
362:   tri->faces       = out.trianglelist;
363:   tri->neighbors   = out.neighborlist;
364:   mesh->holes       = out.holelist;
365:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
366:   PetscLogObjectMemory(mesh, (mesh->numEdges*2+mesh->numEdges+mesh->numFaces*mesh->numCorners+mesh->numFaces*3) * sizeof(int));

368:   return(0);
369: }

371: /*
372:   MeshRefine_Triangle - This function refines a two dimensional unstructured mesh
373:   with Triangle using area constraints.

375:   Collective on Mesh

377:   Input Parameters:
378: + oldMesh - The mesh begin refined
379: - area    - A function which gives an area constraint when evaluated inside an element

381:   Output Parameter:
382: . mesh    - The refined mesh

384:   Level: developer

386: .keywords: mesh, Triangle, refine
387: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
388: */
389: int MeshRefine_Triangle(Mesh oldMesh, PointFunction area, Mesh mesh)
390: {
391:   int                  rank = oldMesh->part->rank;
392:   Mesh_Triangular     *tri;
393:   struct triangulateio in;            /* Input  for Triangle mesh generator */
394:   struct triangulateio out;           /* Output for Triangle mesh generator */
395:   char                 cmd_line[256]; /* Command line for Triangle */
396:   MPI_Comm             comm;
397:   PetscTruth           opt;
398:   int                  ierr;

401:   PetscNew(Mesh_Triangular, &tri);
402:   PetscLogObjectMemory(mesh, sizeof(Mesh_Triangular));
403:   PetscMemcpy(mesh->ops, oldMesh->ops, sizeof(struct _MeshOps));
404:   mesh->data             = (void *) tri;
405:   PetscStrallocpy(MESH_TRIANGULAR_2D,            &mesh->type_name);
406:   PetscStrallocpy(MESH_SER_TRIANGULAR_2D_BINARY, &mesh->serialize_name);
407:   mesh->dim              = 2;
408:   mesh->isPeriodic       = oldMesh->isPeriodic;
409:   mesh->isPeriodicDim[0] = oldMesh->isPeriodicDim[0];
410:   mesh->isPeriodicDim[1] = oldMesh->isPeriodicDim[1];
411:   mesh->isPeriodicDim[2] = oldMesh->isPeriodicDim[2];
412:   mesh->nodeOrdering     = PETSC_NULL;
413:   mesh->partitioned      = 0;
414:   mesh->highlightElement = -1;
415:   mesh->usr              = oldMesh->usr;

417:   /* Copy function list */
418:   if (mesh->qlist != PETSC_NULL) {
419:     PetscFListDestroy(&mesh->qlist);
420:   }
421:   PetscFListDuplicate(oldMesh->qlist, &mesh->qlist);

423:   /* Initialize communication structures for Triangle */
424:   MeshInitRefineInput_Triangle(oldMesh, area, &in);
425:   MeshInitOutput_Triangle(&out);

427:   /* Generate a refined mesh */
428:   if (rank == 0) {
429:     /* Generate the refined mesh and check whether we should create midnodes */
430:     if (oldMesh->numCorners == 3) {
431:       PetscStrcpy(cmd_line, "qcenzQ");
432:     } else if (oldMesh->numCorners == 6) {
433:       PetscStrcpy(cmd_line, "qcenzo2Q");
434:     } else {
435:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", oldMesh->numCorners);
436:     }
437:     PetscOptionsGetString("ref_", "-triangle_cmd", cmd_line, 251, &opt);
438:     PetscStrcat(cmd_line, "pra");
439:     triangulate(cmd_line, &in, &out, PETSC_NULL);

441:     /* Get rid of extra information */
442:     PetscFree(out.segmentlist);
443:     PetscFree(out.segmentmarkerlist);
444:   }
445:   MeshFinalizeRefineInput_Triangle(oldMesh, &in);

447:   /* Communicate values */
448:   MeshDistribute_Triangle(oldMesh, &out);

450:   /* Store the information */
451:   mesh->numBd       = oldMesh->numBd;
452:   mesh->numNodes    = out.numberofpoints;
453:   mesh->numEdges    = out.numberofedges;
454:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
455:   mesh->numFaces    = out.numberoftriangles;
456:   mesh->numCorners  = out.numberofcorners;
457:   tri->nodes       = out.pointlist;
458:   tri->markers     = out.pointmarkerlist;
459:   tri->edges       = out.edgelist;
460:   tri->edgemarkers = out.edgemarkerlist;
461:   tri->faces       = out.trianglelist;
462:   tri->neighbors   = out.neighborlist;
463:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
464:   PetscLogObjectMemory(mesh, (mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners + mesh->numFaces*3) * sizeof(int));
465:   tri->nodesOld    = PETSC_NULL;

467:   /* Copy holes from previous mesh */
468:   mesh->numHoles    = oldMesh->numHoles;
469:   mesh->holes       = PETSC_NULL;
470:   if (mesh->numHoles > 0) {
471:     PetscMalloc(mesh->numHoles*2 * sizeof(double), &mesh->holes);
472:     PetscMemcpy(mesh->holes, oldMesh->holes, mesh->numHoles*2 * sizeof(double));
473:     PetscLogObjectMemory(mesh, oldMesh->numHoles*2 * sizeof(double));
474:   }

476:   /* Calculate maximum degree of vertices */
477:   MeshSetupSupport_Triangular_2D(mesh);

479:   /* Construct derived and boundary information */
480:   MeshSetupBoundary_Triangular_2D(mesh);

482: #ifdef PETSC_USE_BOPT_g
483:   /* Check mesh integrity */
484:   MeshDebug_Triangular_2D(mesh, PETSC_FALSE);
485: #endif

487:   /* Reorder nodes before distributing mesh */
488:   PetscOptionsHasName(oldMesh->prefix, "-mesh_reorder", &opt);
489:   if (opt == PETSC_TRUE) {
490:     /* MUST FIX: Since we have duplicated the whole mesh, we may impersonate a serial mesh */
491:     comm       = mesh->comm;
492:     mesh->comm = PETSC_COMM_SELF;
493:     ierr       = MeshGetOrdering(mesh, MESH_ORDER_TRIANGULAR_2D_RCM, &mesh->nodeOrdering);
494:     mesh->comm = comm;

496:     /* Permute arrays implicitly numbered by node numbers */
497:     AOApplicationToPetscPermuteReal(mesh->nodeOrdering, 2, tri->nodes);
498:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->markers);
499:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->degrees);
500:     /* Renumber arrays dependent on the canonical node numbering */
501:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numEdges*2,                   tri->edges);
502:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numFaces*oldMesh->numCorners, tri->faces);
503:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numBdNodes,                   tri->bdNodes);
504:   }

506:   /* Initialize partition */
507:   MeshPartition(mesh);

509:   return(0);
510: }

512: /*-------------------------------------------- Periodic Interface ---------------------------------------------------*/
513: /* MeshInitInput_Periodic
514:   This function initializes the input structure for Triangle in the periodic case.

516:   Collective on Mesh

518:   Output Parameters:
519: + mesh     - The mesh
520: . inputCtx - The structure for communicating with Triangle
521: - oldHoles - The previous holes, which are saved for MeshFinalizeOutput_Periodic

523:   Level: developer

525: .keywords mesh, periodic, Triangle
526: .seealso MeshFinalizeOutput_Periodic(), MeshInitInput_Triangle()
527: */
528: int MeshInitInput_Periodic(Mesh mesh, struct triangulateio *inputCtx, double **oldHoles)
529: {
530:   int     numPoints = inputCtx->numberofpoints;
531:   double *points    = inputCtx->pointlist;
532:   int     numHoles  = inputCtx->numberofholes;
533:   double *holes     = inputCtx->holelist;
534:   double  r_inner   = mesh->sizeX / (2*PETSC_PI);
535:   double *newHoles;
536:   double  r, theta;
537:   int     rank;
538:   int     p, h;
539:   int     ierr;

542:   MPI_Comm_rank(mesh->comm, &rank);
543:   if (rank == 0) {
544:     /* Map the domain on to an annulus -
545:        We will preserve the periodic length on the inner surface, so that

547:          2 pi r_inner     = mesh->sizeX

549:        and we preserve the bounding rectangle, so that

551:          r_outer - r_inner = mesh->sizeY

553:        where we reverse X and Y if the domain is periodic in Y. Also, we
554:        let theta map along the x-axis in the reverse direction to preserve
555:        the orientation of objects. The equations are therefore

557:          r      = r_inner + (y - mesh->startY)
558:          theta = (2 pi (mesh->sizeX - (x - mesh->startX)) / mesh->sizeX)
559:   */
560:     for(p = 0; p < numPoints; p++) {
561:       r             = r_inner + (points[p*2+1] - mesh->startY);
562:       theta         = 2*PETSC_PI * ((mesh->sizeX - (points[p*2] - mesh->startX)) / mesh->sizeX);
563:       points[p*2]   = r*cos(theta);
564:       points[p*2+1] = r*sin(theta);
565:     }
566:     /* Add the hole in the center of the annulus */
567:     inputCtx->numberofholes = numHoles+1;
568:     PetscMalloc((numHoles+1)*2 * sizeof(double), &newHoles);
569:     for(h = 0; h < numHoles; h++) {
570:       r               = r_inner + (holes[h*2+1] - mesh->startY);
571:       theta           = 2*PETSC_PI * ((mesh->sizeX - (holes[h*2] - mesh->startX)) / mesh->sizeX);
572:       newHoles[h*2]   = r*cos(theta);
573:       newHoles[h*2+1] = r*sin(theta);
574:     }
575:     newHoles[numHoles*2]   = 0.0;
576:     newHoles[numHoles*2+1] = 0.0;
577:     inputCtx->holelist     = newHoles;
578:     if (oldHoles != PETSC_NULL) {
579:       *oldHoles            = holes;
580:     }
581:   }
582:   return(0);
583: }

585: /* MeshFinalizeOutput_Periodic
586:   This function initializes the input structure for Triangle in the periodic case.

588:   Collective on Mesh

590:   Output Parameters:
591: . mesh      - The mesh
592: . outputCtx - The structure for communicating with Triangle
593: . oldHoles  - The previous holes, which were saved in MeshInitInput_Periodic

595:   Level: developer

597: .keywords mesh, periodic, Triangle
598: .seealso MeshInitInput_Periodic()
599: */
600: int MeshFinalizeOutput_Periodic(Mesh mesh, struct triangulateio *outputCtx, double **oldHoles)
601: {
602:   int     numPoints = outputCtx->numberofpoints;
603:   double *points    = outputCtx->pointlist;
604:   int    *markers   = outputCtx->pointmarkerlist;
605:   int     numHoles  = outputCtx->numberofholes;
606:   double  r_inner   = mesh->sizeX / (2*PETSC_PI);
607:   double  r_outer   = r_inner + mesh->sizeY;
608:   double  r, theta;
609:   int     p;
610:   int     ierr;

613:   /* Map the domain from an annulus to a rectangle -
614:        The inverse map is given by

616:          x = mesh->sizeX ((2 pi - theta) / 2 pi) + mesh->startX
617:          y = r - r_inner + mesh->startY
618:   */
619:   for(p = 0; p < numPoints; p++) {
620:     /* If points were inserted on he boundary, we must push them to the correct radius */
621:     if (markers[p] == 2) {        /* TOP_BD    */
622:       r           = r_outer;
623:     } else if (markers[p] == 3) { /* BOTTOM_BD */
624:       r           = r_inner;
625:     } else {
626:       r           = hypot(points[p*2], points[p*2+1]);
627:     }
628:     theta         = atan2(points[p*2+1], points[p*2]);
629:     if (theta < 0)
630:       theta       = 2*PETSC_PI + theta;
631:     points[p*2]   = ((2*PETSC_PI - theta) / (2*PETSC_PI))*mesh->sizeX + mesh->startX;
632:     points[p*2+1] = r - r_inner + mesh->startY;
633:   }
634:   /* Reset the holes */
635:   PetscFree(outputCtx->holelist);
636:   outputCtx->numberofholes = numHoles-1;
637:   if (oldHoles != PETSC_NULL)
638:     outputCtx->holelist    = *oldHoles;
639:   return(0);
640: }

642: /* MeshCreate_Periodic
643:   This function creates a 2D periodic unstructured mesh using Triangle.

645:   Collective on Mesh

647:   Input Parameter:
648: . numCorners  - The number of nodes in an element

650:   Input Parameters from bdCtx:
651: + numBD       - The number of closed boundaries in the geometry, or different markers
652: . numVertices - The umber of boundary points
653: . vertices    - The (x,y) coordinates of the boundary points
654: . markers     - The boundary markers for nodes, 0 indicates an interior point, each boundary must have a different marker
655: . numSegments - The number of boundary segments
656: . segments    - The endpoints of boundary segments or PETSC_NULL
657: . segMarkers  - The boundary markers for each segment
658: . numHoles    - The number of holes
659: - holes       - The (x,y) coordinates of holes or PETSC_NULL

661:   Output Parameter:
662: . mesh        - The new mesh created by Triangle

664:   Level: developer

666: .keywords mesh, Triangle
667: .seealso MeshInitInput_Periodic(), MeshFinalizeOutput_Periodic()
668: */
669: int MeshCreate_Periodic(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
670: {
671:   Mesh_Triangular     *tri = (Mesh_Triangular *) mesh->data;
672:   struct triangulateio in;            /* Input  for Triangle mesh generator */
673:   struct triangulateio out;           /* Output for Triangle mesh generator */
674:   char                 cmd_line[256]; /* Command line for Triangle */
675:   double              *oldHoles;
676:   int                  rank;
677:   PetscTruth           opt;
678:   int                  ierr;

681:   /* Initialize communication structures for Triangle */
682:   MeshInitInput_Triangle(&in);
683:   MeshInitOutput_Triangle(&out);

685:   MPI_Comm_rank(mesh->comm, &rank);
686:   if (rank == 0) {
687:     /* Setup boundaries and holes */
688:     in.numberofpoints      = bdCtx->numVertices;
689:     if (bdCtx->numVertices > 0) {
692:       in.pointlist         = bdCtx->vertices;
693:       in.pointmarkerlist   = bdCtx->markers;
694:     }
695:     in.numberofsegments    = bdCtx->numSegments;
696:     if (bdCtx->numSegments > 0) {
699:       in.segmentlist       = bdCtx->segments;
700:       in.segmentmarkerlist = bdCtx->segMarkers;
701:     }
702:     in.numberofholes       = bdCtx->numHoles;
703:     if (bdCtx->numHoles    > 0) {
704:       /* We keep these */
706:       PetscMalloc(in.numberofholes*2 * sizeof(double), &in.holelist);
707:       PetscMemcpy(in.holelist, bdCtx->holes, in.numberofholes*2 * sizeof(double));
708:     }

710:     /* Generate boundary compatible with Triangle */
711:     MeshInitInput_Periodic(mesh, &in, &oldHoles);

713:     /* Generate the inital coarse mesh and check whether we should create midnodes */
714:     if (numCorners == 3) {
715:       PetscStrcpy(cmd_line, "pqcenzQ");
716:     } else if (numCorners == 6) {
717:       PetscStrcpy(cmd_line, "pqcenzo2Q");
718:     } else {
719:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", numCorners);
720:     }
721:     PetscOptionsGetString(PETSC_NULL, "-triangle_cmd", cmd_line, 255, &opt);
722:     triangulate(cmd_line, &in, &out, PETSC_NULL);

724:     /* Remap output of Triangle to the periodic domain */
725:     MeshFinalizeOutput_Periodic(mesh, &out, &oldHoles);

727:     /* Get rid of extra information */
728:     PetscFree(out.segmentlist);
729:     PetscFree(out.segmentmarkerlist);
730:   }

732:   /* Communicate values */
733:   MeshDistribute_Triangle(mesh, &out);

735:   /* Store the information */
736:   mesh->numBd       = bdCtx->numBd;
737:   mesh->numNodes    = out.numberofpoints;
738:   mesh->numEdges    = out.numberofedges;
739:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
740:   mesh->numFaces    = out.numberoftriangles;
741:   mesh->numCorners  = out.numberofcorners;
742:   mesh->numHoles    = out.numberofholes;
743:   tri->nodes       = out.pointlist;
744:   tri->markers     = out.pointmarkerlist;
745:   tri->edges       = out.edgelist;
746:   tri->edgemarkers = out.edgemarkerlist;
747:   tri->faces       = out.trianglelist;
748:   tri->neighbors   = out.neighborlist;
749:   mesh->holes       = out.holelist;
750:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
751:   PetscLogObjectMemory(mesh, (mesh->numEdges*2+mesh->numEdges+mesh->numFaces*mesh->numCorners+mesh->numFaces*3) * sizeof(int));

753:   return(0);
754: }

756: /*
757:   MeshRefine_Periodic - This function refines a two dimensional periodic unstructured mesh
758:   with Triangle using area constraints.

760:   Collective on Mesh

762:   Input Parameters:
763: + oldMesh - The mesh begin refined
764: - area    - A function which gives an area constraint when evaluated inside an element

766:   Output Parameter:
767: . mesh    - The refined mesh

769:   Level: developer

771: .keywords: mesh, Triangle, refine, periodic
772: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
773: */
774: int MeshRefine_Periodic(Mesh oldMesh, PointFunction area, Mesh mesh)
775: {
776:   int                  rank = oldMesh->part->rank;
777:   Mesh_Triangular     *tri;
778:   struct triangulateio in;            /* Input  for Triangle mesh generator */
779:   struct triangulateio out;           /* Output for Triangle mesh generator */
780:   char                 cmd_line[256]; /* Command line for Triangle */
781:   MPI_Comm             comm;
782:   PetscTruth           opt;
783:   int                  ierr;

786:   PetscNew(Mesh_Triangular, &tri);
787:   PetscLogObjectMemory(mesh, sizeof(Mesh_Triangular));
788:   PetscMemcpy(mesh->ops, oldMesh->ops, sizeof(struct _MeshOps));
789:   mesh->data             = (void *) tri;
790:   PetscStrallocpy(MESH_TRIANGULAR_2D,            &mesh->type_name);
791:   PetscStrallocpy(MESH_SER_TRIANGULAR_2D_BINARY, &mesh->serialize_name);
792:   mesh->dim              = 2;
793:   mesh->isPeriodic       = oldMesh->isPeriodic;
794:   mesh->isPeriodicDim[0] = oldMesh->isPeriodicDim[0];
795:   mesh->isPeriodicDim[1] = oldMesh->isPeriodicDim[1];
796:   mesh->isPeriodicDim[2] = oldMesh->isPeriodicDim[2];
797:   mesh->nodeOrdering     = PETSC_NULL;
798:   mesh->partitioned      = 0;
799:   mesh->highlightElement = -1;
800:   mesh->usr              = oldMesh->usr;

802:   /* Copy function list */
803:   if (mesh->qlist != PETSC_NULL) {
804:     PetscFListDestroy(&mesh->qlist);
805:   }
806:   PetscFListDuplicate(oldMesh->qlist, &mesh->qlist);

808:   /* Initialize communication structures for Triangle */
809:   MeshInitRefineInput_Triangle(oldMesh, area, &in);
810:   MeshInitInput_Periodic(oldMesh, &in, PETSC_NULL);
811:   MeshInitOutput_Triangle(&out);

813:   /* Copy extent from previous mesh: don't make separate boxes for now */
814:   mesh->startX    = oldMesh->startX;
815:   mesh->endX      = oldMesh->endX;
816:   mesh->startY    = oldMesh->startY;
817:   mesh->endY      = oldMesh->endY;
818:   mesh->locStartX = mesh->startX;
819:   mesh->locEndX   = mesh->endX;
820:   mesh->locStartY = mesh->startY;
821:   mesh->locEndY   = mesh->endY;
822:   mesh->sizeX     = mesh->endX    - mesh->startX;
823:   mesh->sizeY     = mesh->endY    - mesh->startY;
824:   mesh->locSizeX  = mesh->locEndX - mesh->locStartX;
825:   mesh->locSizeY  = mesh->locEndY - mesh->locStartY;

827:   /* Generate a refined mesh */
828:   if (rank == 0) {
829:     /* Generate the refined mesh and check whether we should create midnodes */
830:     if (oldMesh->numCorners == 3) {
831:       PetscStrcpy(cmd_line, "qcenzQ");
832:     } else if (oldMesh->numCorners == 6) {
833:       PetscStrcpy(cmd_line, "qcenzo2Q");
834:     } else {
835:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", oldMesh->numCorners);
836:     }
837:     PetscOptionsGetString("ref_", "-triangle_cmd", cmd_line, 251, &opt);
838:     PetscStrcat(cmd_line, "pra");
839:     triangulate(cmd_line, &in, &out, PETSC_NULL);

841:     /* Remap output of Triangle to the periodic domain */
842:     MeshFinalizeOutput_Periodic(mesh, &out, PETSC_NULL);

844:     /* Get rid of extra information */
845:     PetscFree(out.segmentlist);
846:     PetscFree(out.segmentmarkerlist);
847:   }
848:   MeshFinalizeRefineInput_Triangle(oldMesh, &in);

850:   /* Communicate values */
851:   MeshDistribute_Triangle(oldMesh, &out);

853:   /* Store the information */
854:   mesh->numBd       = oldMesh->numBd;
855:   mesh->numNodes    = out.numberofpoints;
856:   mesh->numEdges    = out.numberofedges;
857:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
858:   mesh->numFaces    = out.numberoftriangles;
859:   mesh->numCorners  = out.numberofcorners;
860:   tri->nodes       = out.pointlist;
861:   tri->markers     = out.pointmarkerlist;
862:   tri->edges       = out.edgelist;
863:   tri->edgemarkers = out.edgemarkerlist;
864:   tri->faces       = out.trianglelist;
865:   tri->neighbors   = out.neighborlist;
866:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
867:   PetscLogObjectMemory(mesh, (mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners + mesh->numFaces*3) * sizeof(int));
868:   tri->nodesOld    = PETSC_NULL;

870:   /* Copy holes from previous mesh */
871:   mesh->numHoles    = oldMesh->numHoles;
872:   mesh->holes       = PETSC_NULL;
873:   if (mesh->numHoles > 0) {
874:     PetscMalloc(mesh->numHoles*2 * sizeof(double), &mesh->holes);
875:     PetscMemcpy(mesh->holes, oldMesh->holes, mesh->numHoles*2 * sizeof(double));
876:     PetscLogObjectMemory(mesh, oldMesh->numHoles*2 * sizeof(double));
877:   }

879:   /* Calculate maximum degree of vertices */
880:   MeshSetupSupport_Triangular_2D(mesh);

882:   /* Construct derived and boundary information */
883:   MeshSetupBoundary_Triangular_2D(mesh);

885: #ifdef PETSC_USE_BOPT_g
886:   /* Check mesh integrity */
887:   MeshDebug_Triangular_2D(mesh, PETSC_FALSE);
888: #endif

890:   /* Reorder nodes before distributing mesh */
891:   PetscOptionsHasName(oldMesh->prefix, "-mesh_reorder", &opt);
892:   if (opt == PETSC_TRUE) {
893:     /* MUST FIX: Since we have duplicated the whole mesh, we may impersonate a serial mesh */
894:     comm       = mesh->comm;
895:     mesh->comm = PETSC_COMM_SELF;
896:     ierr       = MeshGetOrdering(mesh, MESH_ORDER_TRIANGULAR_2D_RCM, &mesh->nodeOrdering);
897:     mesh->comm = comm;

899:     /* Permute arrays implicitly numbered by node numbers */
900:     AOApplicationToPetscPermuteReal(mesh->nodeOrdering, 2, tri->nodes);
901:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->markers);
902:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->degrees);
903:     /* Renumber arrays dependent on the canonical node numbering */
904:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numEdges*2,                   tri->edges);
905:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numFaces*oldMesh->numCorners, tri->faces);
906:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numBdNodes,                   tri->bdNodes);
907:   }

909:   /* Initialize partition */
910:   MeshPartition(mesh);

912:   /* Reset the midnodes */
913:   MeshResetNodes(mesh, PETSC_TRUE);
914:   return(0);
915: }

917: #endif /* PETSC_HAVE_TRIANGLE */