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: #undef  __FUNCT__
 21: /* MeshInitInput_Triangle
 22:   This function initializes the input structure for Triangle.

 24:   Collective on Mesh

 26:   Output Parameter:
 27: . inputCtx - Structure for communicating with Triangle

 29:   Level: developer

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

 53: #undef  __FUNCT__
 55: /* MeshInitRefineInput_Triangle
 56:   This function fills the input structure for Triangle with mesh information before refinement.

 58:   Collective on Mesh

 60:   Input Parameters:
 61: + mesh    - The mesh being refined
 62: - area     - A function which gives an area constraint when evaluated inside an element

 64:   Output Parameter:
 65: . inputCtx - Structure for communicating with Triangle

 67:   Level: developer

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

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

 91:   if (rank == 0) {
 92:     PetscMalloc(mesh->numBdEdges*2 * sizeof(int), &segments);
 93:     PetscMalloc(mesh->numBdEdges   * sizeof(int), &segmarkers);

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

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

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

141:     /* Create list of area constraints */
142:     PetscMalloc(numFaces * sizeof(double), &inputCtx->trianglearealist);
143:     for(elem = 0; elem < numFaces; elem++) {
144:       PetscReal maxArea;

146:       /* Calculate the centroid of the element */
147:       x = MeshPeriodicX(mesh, (nodes[faces[elem*numCorners]*2]   + nodes[faces[elem*numCorners+1]*2]   +
148:                         nodes[faces[elem*numCorners+2]*2])  /3.0);
149:       y = MeshPeriodicY(mesh, (nodes[faces[elem*numCorners]*2+1] + nodes[faces[elem*numCorners+1]*2+1] +
150:                         nodes[faces[elem*numCorners+2]*2+1])/3.0);
151:       (*area)(1, 1, &x, &y, PETSC_NULL, &maxArea, mesh->usr);
152:       inputCtx->trianglearealist[elem] = maxArea;
153:     }
154:   }
155: #ifdef PETSC_USE_BOPT_g
156:   PetscTrValid(__LINE__, __FUNCT__, __FILE__, __SDIR__);
157: #endif

159:   return(0);
160: }

162: #undef  __FUNCT__
164: /* MeshFinalizeRefineInput_Triangle
165:   This function cleans up the input structure for Triangle after refinement.

167:   Collective on Mesh

169:   Input Parameters:
170: . mesh     - The mesh being refined

172:   Output Parameters:
173: . inputCtx - Structure for communicating with Triangle

175:   Level: developer

177: .keywords mesh, Triangle
178: .seealso MeshInitRefineInput_Triangle(), MeshInitInput_Triangle()
179: */
180: int MeshFinalizeRefineInput_Triangle(Mesh mesh, struct triangulateio *inputCtx)
181: {

185:   if (mesh->part->rank == 0) {
186:     PetscFree(inputCtx->pointlist);
187:     PetscFree(inputCtx->pointmarkerlist);
188:     PetscFree(inputCtx->trianglelist);
189:     PetscFree(inputCtx->segmentlist);
190:     PetscFree(inputCtx->segmentmarkerlist);
191:     PetscFree(inputCtx->trianglearealist);
192:   }
193:   return(0);
194: }

196: #undef  __FUNCT__
198: /* MeshInitOutput_Triangle
199:   This function initializes the output structure for Triangle.

201:   Collective on Mesh

203:   Output Parameter:
204: . outputCtx - Structure for communicating with Triangle

206:   Level: developer

208: .keywords mesh, Triangle
209: .seealso MeshInitInput_Triangle()
210: */
211: int MeshInitOutput_Triangle(struct triangulateio *outputCtx)
212: {
214:   /* Setup the output structure */
215:   outputCtx->pointlist             = PETSC_NULL;
216:   outputCtx->pointattributelist    = PETSC_NULL;
217:   outputCtx->pointmarkerlist       = PETSC_NULL;
218:   outputCtx->trianglelist          = PETSC_NULL;
219:   outputCtx->triangleattributelist = PETSC_NULL;
220:   outputCtx->neighborlist          = PETSC_NULL;
221:   outputCtx->segmentlist           = PETSC_NULL;
222:   outputCtx->segmentmarkerlist     = PETSC_NULL;
223:   outputCtx->edgelist              = PETSC_NULL;
224:   outputCtx->edgemarkerlist        = PETSC_NULL;
225:   return(0);
226: }

228: #undef  __FUNCT__
230: /* MeshDistribute_Triangle
231:   This function distributes the output from Triangle identically among processors.

233:   Collective on Mesh

235:   Output Parameters:
236: + mesh - The mesh
237: - out  - Structure for communicating with Triangle

239:   Level: developer

241: .keywords mesh, Triangle
242: .seealso MeshInitInput_Triangle
243: */
244: int MeshDistribute_Triangle(Mesh mesh, struct triangulateio *out)
245: {
246:   int rank;

250:   /* Communicate values */
251:   MPI_Comm_rank(mesh->comm, &rank);
252:   MPI_Bcast(&out->numberofpoints,    1, MPI_INT, 0, mesh->comm);
253:   MPI_Bcast(&out->numberofedges,     1, MPI_INT, 0, mesh->comm);
254:   MPI_Bcast(&out->numberoftriangles, 1, MPI_INT, 0, mesh->comm);
255:   MPI_Bcast(&out->numberofcorners,   1, MPI_INT, 0, mesh->comm);
256:   if (rank != 0) {
257:     out->numberofholes   = 0;
258:     out->holelist        = PETSC_NULL;
259:     PetscMalloc(out->numberofpoints*2                       * sizeof(double), &out->pointlist);
260:     PetscMalloc(out->numberofpoints                         * sizeof(int),    &out->pointmarkerlist);
261:     PetscMalloc(out->numberofedges*2                        * sizeof(int),    &out->edgelist);
262:     PetscMalloc(out->numberofedges                          * sizeof(int),    &out->edgemarkerlist);
263:     PetscMalloc(out->numberoftriangles*3                    * sizeof(int),    &out->neighborlist);
264:     PetscMalloc(out->numberoftriangles*out->numberofcorners * sizeof(int),    &out->trianglelist);
265:   }
266:   MPI_Bcast(out->pointlist,       out->numberofpoints*2,    MPI_DOUBLE, 0, mesh->comm);
267:   MPI_Bcast(out->pointmarkerlist, out->numberofpoints,      MPI_INT,    0, mesh->comm);
268:   MPI_Bcast(out->edgelist,        out->numberofedges*2,     MPI_INT,    0, mesh->comm);
269:   MPI_Bcast(out->edgemarkerlist,  out->numberofedges,       MPI_INT,    0, mesh->comm);
270:   MPI_Bcast(out->neighborlist,    out->numberoftriangles*3, MPI_INT,    0, mesh->comm);
271:   MPI_Bcast(out->trianglelist,    out->numberoftriangles*out->numberofcorners, MPI_INT, 0, mesh->comm);
272:   return(0);
273: }

275: #undef  __FUNCT__
277: /* MeshCreate_Triangle
278:   This function creates a 2D unstructured mesh using Triangle.

280:   Collective on Mesh

282:   Input Parameter:
283: . numCorners  - The number of nodes in an element

285:   Input Parameters from bdCtx:
286: + numBD       - The number of closed boundaries in the geometry, or different markers
287: . numVertices - The number of boundary points
288: . vertices    - The (x,y) coordinates of the boundary points
289: . markers     - The boundary markers for nodes, 0 indicates an interior point, each boundary must have a different marker
290: . numSegments - The number of boundary segments
291: . segments    - The endpoints of boundary segments or PETSC_NULL
292: . segMarkers  - The boundary markers for each segment
293: . numHoles    - The number of holes
294: - holes       - The (x,y) coordinates of holes or PETSC_NULL

296:   Output Parameter:
297: . mesh        - The new mesh created by Triangle

299:   Level: developer

301: .keywords mesh, Triangle
302: .seealso MeshInitInput_Triangle
303: */
304: int MeshCreate_Triangle(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
305: {
306:   Mesh_Triangular     *tri = (Mesh_Triangular *) mesh->data;
307:   struct triangulateio in;            /* Input  for Triangle mesh generator */
308:   struct triangulateio out;           /* Output for Triangle mesh generator */
309:   char                 cmd_line[256]; /* Command line for Triangle */
310:   int                  rank;
311:   PetscTruth           opt;
312:   int                  ierr;

315:   /* Initialize communication structures for Triangle */
316:   MeshInitInput_Triangle(&in);
317:   MeshInitOutput_Triangle(&out);

319:   MPI_Comm_rank(mesh->comm, &rank);
320:   if (rank == 0) {
321:     /* Setup boundaries and holes */
322:     in.numberofpoints      = bdCtx->numVertices;
323:     if (bdCtx->numVertices > 0) {
326:       in.pointlist         = bdCtx->vertices;
327:       in.pointmarkerlist   = bdCtx->markers;
328:     }
329:     in.numberofsegments    = bdCtx->numSegments;
330:     if (bdCtx->numSegments > 0) {
333:       in.segmentlist       = bdCtx->segments;
334:       in.segmentmarkerlist = bdCtx->segMarkers;
335:     }
336:     in.numberofholes       = bdCtx->numHoles;
337:     if (bdCtx->numHoles    > 0) {
338:       /* We keep these */
340:       PetscMalloc(in.numberofholes*2 * sizeof(double), &in.holelist);
341:       PetscMemcpy(in.holelist, bdCtx->holes, in.numberofholes*2 * sizeof(double));
342:     }

344:     /* Generate the inital coarse mesh and check whether we should create midnodes */
345:     if (numCorners == 3) {
346:       PetscStrcpy(cmd_line, "pqcenzQ");
347:     } else if (numCorners == 6) {
348:       PetscStrcpy(cmd_line, "pqcenzo2Q");
349:     } else {
350:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", numCorners);
351:     }
352:     PetscOptionsGetString(PETSC_NULL, "-triangle_cmd", cmd_line, 255, &opt);
353:     triangulate(cmd_line, &in, &out, PETSC_NULL);

355:     /* Get rid of extra information */
356:     PetscFree(out.segmentlist);
357:     PetscFree(out.segmentmarkerlist);
358:   }

360:   /* Communicate values */
361:   MeshDistribute_Triangle(mesh, &out);

363:   /* Store the information */
364:   mesh->numBd       = bdCtx->numBd;
365:   mesh->numNodes    = out.numberofpoints;
366:   mesh->numEdges    = out.numberofedges;
367:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
368:   mesh->numFaces    = out.numberoftriangles;
369:   mesh->numCorners  = out.numberofcorners;
370:   mesh->numHoles    = out.numberofholes;
371:   tri->nodes       = out.pointlist;
372:   tri->markers     = out.pointmarkerlist;
373:   tri->edges       = out.edgelist;
374:   tri->edgemarkers = out.edgemarkerlist;
375:   tri->faces       = out.trianglelist;
376:   tri->neighbors   = out.neighborlist;
377:   mesh->holes       = out.holelist;
378:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
379:   PetscLogObjectMemory(mesh, (mesh->numEdges*2+mesh->numEdges+mesh->numFaces*mesh->numCorners+mesh->numFaces*3) * sizeof(int));

381:   return(0);
382: }

384: #undef  __FUNCT__
386: /*
387:   MeshRefine_Triangle - This function refines a two dimensional unstructured mesh
388:   with Triangle using area constraints.

390:   Collective on Mesh

392:   Input Parameters:
393: + oldMesh - The mesh begin refined
394: - area    - A function which gives an area constraint when evaluated inside an element

396:   Output Parameter:
397: . mesh    - The refined mesh

399:   Level: developer

401: .keywords: mesh, Triangle, refine
402: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
403: */
404: int MeshRefine_Triangle(Mesh oldMesh, PointFunction area, Mesh mesh)
405: {
406:   int                  rank = oldMesh->part->rank;
407:   Mesh_Triangular     *tri;
408:   struct triangulateio in;            /* Input  for Triangle mesh generator */
409:   struct triangulateio out;           /* Output for Triangle mesh generator */
410:   char                 cmd_line[256]; /* Command line for Triangle */
411:   MPI_Comm             comm;
412:   PetscTruth           opt;
413:   int                  ierr;

416:   PetscNew(Mesh_Triangular, &tri);
417:   PetscLogObjectMemory(mesh, sizeof(Mesh_Triangular));
418:   PetscMemcpy(mesh->ops, oldMesh->ops, sizeof(struct _MeshOps));
419:   mesh->data             = (void *) tri;
420:   PetscStrallocpy(MESH_TRIANGULAR_2D,            &mesh->type_name);
421:   PetscStrallocpy(MESH_SER_TRIANGULAR_2D_BINARY, &mesh->serialize_name);
422:   mesh->dim              = 2;
423:   mesh->isPeriodic       = oldMesh->isPeriodic;
424:   mesh->isPeriodicDim[0] = oldMesh->isPeriodicDim[0];
425:   mesh->isPeriodicDim[1] = oldMesh->isPeriodicDim[1];
426:   mesh->isPeriodicDim[2] = oldMesh->isPeriodicDim[2];
427:   mesh->nodeOrdering     = PETSC_NULL;
428:   mesh->partitioned      = 0;
429:   mesh->highlightElement = -1;
430:   mesh->usr              = oldMesh->usr;

432:   /* Copy function list */
433:   if (mesh->qlist != PETSC_NULL) {
434:     PetscFListDestroy(&mesh->qlist);
435:   }
436:   PetscFListDuplicate(oldMesh->qlist, &mesh->qlist);

438:   /* Initialize communication structures for Triangle */
439:   MeshInitRefineInput_Triangle(oldMesh, area, &in);
440:   MeshInitOutput_Triangle(&out);

442:   /* Generate a refined mesh */
443:   if (rank == 0) {
444:     /* Generate the refined mesh and check whether we should create midnodes */
445:     if (oldMesh->numCorners == 3) {
446:       PetscStrcpy(cmd_line, "qcenzQ");
447:     } else if (oldMesh->numCorners == 6) {
448:       PetscStrcpy(cmd_line, "qcenzo2Q");
449:     } else {
450:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", oldMesh->numCorners);
451:     }
452:     PetscOptionsGetString("ref_", "-triangle_cmd", cmd_line, 251, &opt);
453:     PetscStrcat(cmd_line, "pra");
454:     triangulate(cmd_line, &in, &out, PETSC_NULL);

456:     /* Get rid of extra information */
457:     PetscFree(out.segmentlist);
458:     PetscFree(out.segmentmarkerlist);
459:   }
460:   MeshFinalizeRefineInput_Triangle(oldMesh, &in);

462:   /* Communicate values */
463:   MeshDistribute_Triangle(oldMesh, &out);

465:   /* Store the information */
466:   mesh->numBd       = oldMesh->numBd;
467:   mesh->numNodes    = out.numberofpoints;
468:   mesh->numEdges    = out.numberofedges;
469:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
470:   mesh->numFaces    = out.numberoftriangles;
471:   mesh->numCorners  = out.numberofcorners;
472:   tri->nodes       = out.pointlist;
473:   tri->markers     = out.pointmarkerlist;
474:   tri->edges       = out.edgelist;
475:   tri->edgemarkers = out.edgemarkerlist;
476:   tri->faces       = out.trianglelist;
477:   tri->neighbors   = out.neighborlist;
478:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
479:   PetscLogObjectMemory(mesh, (mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners + mesh->numFaces*3) * sizeof(int));
480:   tri->nodesOld    = PETSC_NULL;

482:   /* Copy holes from previous mesh */
483:   mesh->numHoles    = oldMesh->numHoles;
484:   mesh->holes       = PETSC_NULL;
485:   if (mesh->numHoles > 0) {
486:     PetscMalloc(mesh->numHoles*2 * sizeof(double), &mesh->holes);
487:     PetscMemcpy(mesh->holes, oldMesh->holes, mesh->numHoles*2 * sizeof(double));
488:     PetscLogObjectMemory(mesh, oldMesh->numHoles*2 * sizeof(double));
489:   }

491:   /* Calculate maximum degree of vertices */
492:   MeshSetupSupport_Triangular_2D(mesh);

494:   /* Construct derived and boundary information */
495:   MeshSetupBoundary_Triangular_2D(mesh);

497: #ifdef PETSC_USE_BOPT_g
498:   /* Check mesh integrity */
499:   MeshDebug_Triangular_2D(mesh, PETSC_FALSE);
500: #endif

502:   /* Reorder nodes before distributing mesh */
503:   PetscOptionsHasName(oldMesh->prefix, "-mesh_reorder", &opt);
504:   if (opt == PETSC_TRUE) {
505:     /* MUST FIX: Since we have duplicated the whole mesh, we may impersonate a serial mesh */
506:     comm       = mesh->comm;
507:     mesh->comm = PETSC_COMM_SELF;
508:     ierr       = MeshGetOrdering(mesh, MESH_ORDER_TRIANGULAR_2D_RCM, &mesh->nodeOrdering);
509:     mesh->comm = comm;

511:     /* Permute arrays implicitly numbered by node numbers */
512:     AOApplicationToPetscPermuteReal(mesh->nodeOrdering, 2, tri->nodes);
513:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->markers);
514:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->degrees);
515:     /* Renumber arrays dependent on the canonical node numbering */
516:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numEdges*2,                   tri->edges);
517:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numFaces*oldMesh->numCorners, tri->faces);
518:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numBdNodes,                   tri->bdNodes);
519:   }

521:   /* Initialize partition */
522:   MeshPartition(mesh);

524:   return(0);
525: }

527: /*-------------------------------------------- Periodic Interface ---------------------------------------------------*/
528: #undef  __FUNCT__
530: /* MeshInitInput_Periodic
531:   This function initializes the input structure for Triangle in the periodic case.

533:   Collective on Mesh

535:   Output Parameters:
536: + mesh     - The mesh
537: . inputCtx - The structure for communicating with Triangle
538: - oldHoles - The previous holes, which are saved for MeshFinalizeOutput_Periodic

540:   Level: developer

542: .keywords mesh, periodic, Triangle
543: .seealso MeshFinalizeOutput_Periodic(), MeshInitInput_Triangle()
544: */
545: int MeshInitInput_Periodic(Mesh mesh, struct triangulateio *inputCtx, double **oldHoles)
546: {
547:   int     numPoints = inputCtx->numberofpoints;
548:   double *points    = inputCtx->pointlist;
549:   int     numHoles  = inputCtx->numberofholes;
550:   double *holes     = inputCtx->holelist;
551:   double  r_inner   = mesh->sizeX / (2*PETSC_PI);
552:   double *newHoles;
553:   double  r, theta;
554:   int     rank;
555:   int     p, h;
556:   int     ierr;

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

564:          2 pi r_inner     = mesh->sizeX

566:        and we preserve the bounding rectangle, so that

568:          r_outer - r_inner = mesh->sizeY

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

574:          r      = r_inner + (y - mesh->startY)
575:          theta = (2 pi (mesh->sizeX - (x - mesh->startX)) / mesh->sizeX)
576:   */
577:     for(p = 0; p < numPoints; p++) {
578:       r             = r_inner + (points[p*2+1] - mesh->startY);
579:       theta         = 2*PETSC_PI * ((mesh->sizeX - (points[p*2] - mesh->startX)) / mesh->sizeX);
580:       points[p*2]   = r*cos(theta);
581:       points[p*2+1] = r*sin(theta);
582:     }
583:     /* Add the hole in the center of the annulus */
584:     inputCtx->numberofholes = numHoles+1;
585:     PetscMalloc((numHoles+1)*2 * sizeof(double), &newHoles);
586:     for(h = 0; h < numHoles; h++) {
587:       r               = r_inner + (holes[h*2+1] - mesh->startY);
588:       theta           = 2*PETSC_PI * ((mesh->sizeX - (holes[h*2] - mesh->startX)) / mesh->sizeX);
589:       newHoles[h*2]   = r*cos(theta);
590:       newHoles[h*2+1] = r*sin(theta);
591:     }
592:     newHoles[numHoles*2]   = 0.0;
593:     newHoles[numHoles*2+1] = 0.0;
594:     inputCtx->holelist     = newHoles;
595:     if (oldHoles != PETSC_NULL) {
596:       *oldHoles            = holes;
597:     }
598:   }
599:   return(0);
600: }

602: #undef  __FUNCT__
604: /* MeshFinalizeOutput_Periodic
605:   This function initializes the input structure for Triangle in the periodic case.

607:   Collective on Mesh

609:   Output Parameters:
610: . mesh      - The mesh
611: . outputCtx - The structure for communicating with Triangle
612: . oldHoles  - The previous holes, which were saved in MeshInitInput_Periodic

614:   Level: developer

616: .keywords mesh, periodic, Triangle
617: .seealso MeshInitInput_Periodic()
618: */
619: int MeshFinalizeOutput_Periodic(Mesh mesh, struct triangulateio *outputCtx, double **oldHoles)
620: {
621:   int     numPoints = outputCtx->numberofpoints;
622:   double *points    = outputCtx->pointlist;
623:   int    *markers   = outputCtx->pointmarkerlist;
624:   int     numHoles  = outputCtx->numberofholes;
625:   double  r_inner   = mesh->sizeX / (2*PETSC_PI);
626:   double  r_outer   = r_inner + mesh->sizeY;
627:   double  r, theta;
628:   int     p;
629:   int     ierr;

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

635:          x = mesh->sizeX ((2 pi - theta) / 2 pi) + mesh->startX
636:          y = r - r_inner + mesh->startY
637:   */
638:   for(p = 0; p < numPoints; p++) {
639:     /* If points were inserted on he boundary, we must push them to the correct radius */
640:     if (markers[p] == 2) {        /* TOP_BD    */
641:       r           = r_outer;
642:     } else if (markers[p] == 3) { /* BOTTOM_BD */
643:       r           = r_inner;
644:     } else {
645:       r           = hypot(points[p*2], points[p*2+1]);
646:     }
647:     theta         = atan2(points[p*2+1], points[p*2]);
648:     if (theta < 0)
649:       theta       = 2*PETSC_PI + theta;
650:     points[p*2]   = ((2*PETSC_PI - theta) / (2*PETSC_PI))*mesh->sizeX + mesh->startX;
651:     points[p*2+1] = r - r_inner + mesh->startY;
652:   }
653:   /* Reset the holes */
654:   PetscFree(outputCtx->holelist);
655:   outputCtx->numberofholes = numHoles-1;
656:   if (oldHoles != PETSC_NULL)
657:     outputCtx->holelist    = *oldHoles;
658:   return(0);
659: }

661: #undef  __FUNCT__
663: /* MeshCreate_Periodic
664:   This function creates a 2D periodic unstructured mesh using Triangle.

666:   Collective on Mesh

668:   Input Parameter:
669: . numCorners  - The number of nodes in an element

671:   Input Parameters from bdCtx:
672: + numBD       - The number of closed boundaries in the geometry, or different markers
673: . numVertices - The umber of boundary points
674: . vertices    - The (x,y) coordinates of the boundary points
675: . markers     - The boundary markers for nodes, 0 indicates an interior point, each boundary must have a different marker
676: . numSegments - The number of boundary segments
677: . segments    - The endpoints of boundary segments or PETSC_NULL
678: . segMarkers  - The boundary markers for each segment
679: . numHoles    - The number of holes
680: - holes       - The (x,y) coordinates of holes or PETSC_NULL

682:   Output Parameter:
683: . mesh        - The new mesh created by Triangle

685:   Level: developer

687: .keywords mesh, Triangle
688: .seealso MeshInitInput_Periodic(), MeshFinalizeOutput_Periodic()
689: */
690: int MeshCreate_Periodic(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
691: {
692:   Mesh_Triangular     *tri = (Mesh_Triangular *) mesh->data;
693:   struct triangulateio in;            /* Input  for Triangle mesh generator */
694:   struct triangulateio out;           /* Output for Triangle mesh generator */
695:   char                 cmd_line[256]; /* Command line for Triangle */
696:   double              *oldHoles;
697:   int                  rank;
698:   PetscTruth           opt;
699:   int                  ierr;

702:   /* Initialize communication structures for Triangle */
703:   MeshInitInput_Triangle(&in);
704:   MeshInitOutput_Triangle(&out);

706:   MPI_Comm_rank(mesh->comm, &rank);
707:   if (rank == 0) {
708:     /* Setup boundaries and holes */
709:     in.numberofpoints      = bdCtx->numVertices;
710:     if (bdCtx->numVertices > 0) {
713:       in.pointlist         = bdCtx->vertices;
714:       in.pointmarkerlist   = bdCtx->markers;
715:     }
716:     in.numberofsegments    = bdCtx->numSegments;
717:     if (bdCtx->numSegments > 0) {
720:       in.segmentlist       = bdCtx->segments;
721:       in.segmentmarkerlist = bdCtx->segMarkers;
722:     }
723:     in.numberofholes       = bdCtx->numHoles;
724:     if (bdCtx->numHoles    > 0) {
725:       /* We keep these */
727:       PetscMalloc(in.numberofholes*2 * sizeof(double), &in.holelist);
728:       PetscMemcpy(in.holelist, bdCtx->holes, in.numberofholes*2 * sizeof(double));
729:     }

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

734:     /* Generate the inital coarse mesh and check whether we should create midnodes */
735:     if (numCorners == 3) {
736:       PetscStrcpy(cmd_line, "pqcenzQ");
737:     } else if (numCorners == 6) {
738:       PetscStrcpy(cmd_line, "pqcenzo2Q");
739:     } else {
740:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", numCorners);
741:     }
742:     PetscOptionsGetString(PETSC_NULL, "-triangle_cmd", cmd_line, 255, &opt);
743:     triangulate(cmd_line, &in, &out, PETSC_NULL);

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

748:     /* Get rid of extra information */
749:     PetscFree(out.segmentlist);
750:     PetscFree(out.segmentmarkerlist);
751:   }

753:   /* Communicate values */
754:   MeshDistribute_Triangle(mesh, &out);

756:   /* Store the information */
757:   mesh->numBd       = bdCtx->numBd;
758:   mesh->numNodes    = out.numberofpoints;
759:   mesh->numEdges    = out.numberofedges;
760:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
761:   mesh->numFaces    = out.numberoftriangles;
762:   mesh->numCorners  = out.numberofcorners;
763:   mesh->numHoles    = out.numberofholes;
764:   tri->nodes       = out.pointlist;
765:   tri->markers     = out.pointmarkerlist;
766:   tri->edges       = out.edgelist;
767:   tri->edgemarkers = out.edgemarkerlist;
768:   tri->faces       = out.trianglelist;
769:   tri->neighbors   = out.neighborlist;
770:   mesh->holes       = out.holelist;
771:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
772:   PetscLogObjectMemory(mesh, (mesh->numEdges*2+mesh->numEdges+mesh->numFaces*mesh->numCorners+mesh->numFaces*3) * sizeof(int));

774:   return(0);
775: }

777: #undef  __FUNCT__
779: /*
780:   MeshRefine_Periodic - This function refines a two dimensional periodic unstructured mesh
781:   with Triangle using area constraints.

783:   Collective on Mesh

785:   Input Parameters:
786: + oldMesh - The mesh begin refined
787: - area    - A function which gives an area constraint when evaluated inside an element

789:   Output Parameter:
790: . mesh    - The refined mesh

792:   Level: developer

794: .keywords: mesh, Triangle, refine, periodic
795: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
796: */
797: int MeshRefine_Periodic(Mesh oldMesh, PointFunction area, Mesh mesh)
798: {
799:   int                  rank = oldMesh->part->rank;
800:   Mesh_Triangular     *tri;
801:   struct triangulateio in;            /* Input  for Triangle mesh generator */
802:   struct triangulateio out;           /* Output for Triangle mesh generator */
803:   char                 cmd_line[256]; /* Command line for Triangle */
804:   MPI_Comm             comm;
805:   PetscTruth           opt;
806:   int                  ierr;

809:   PetscNew(Mesh_Triangular, &tri);
810:   PetscLogObjectMemory(mesh, sizeof(Mesh_Triangular));
811:   PetscMemcpy(mesh->ops, oldMesh->ops, sizeof(struct _MeshOps));
812:   mesh->data             = (void *) tri;
813:   PetscStrallocpy(MESH_TRIANGULAR_2D,            &mesh->type_name);
814:   PetscStrallocpy(MESH_SER_TRIANGULAR_2D_BINARY, &mesh->serialize_name);
815:   mesh->dim              = 2;
816:   mesh->isPeriodic       = oldMesh->isPeriodic;
817:   mesh->isPeriodicDim[0] = oldMesh->isPeriodicDim[0];
818:   mesh->isPeriodicDim[1] = oldMesh->isPeriodicDim[1];
819:   mesh->isPeriodicDim[2] = oldMesh->isPeriodicDim[2];
820:   mesh->nodeOrdering     = PETSC_NULL;
821:   mesh->partitioned      = 0;
822:   mesh->highlightElement = -1;
823:   mesh->usr              = oldMesh->usr;

825:   /* Copy function list */
826:   if (mesh->qlist != PETSC_NULL) {
827:     PetscFListDestroy(&mesh->qlist);
828:   }
829:   PetscFListDuplicate(oldMesh->qlist, &mesh->qlist);

831:   /* Initialize communication structures for Triangle */
832:   MeshInitRefineInput_Triangle(oldMesh, area, &in);
833:   MeshInitInput_Periodic(oldMesh, &in, PETSC_NULL);
834:   MeshInitOutput_Triangle(&out);

836:   /* Copy extent from previous mesh: don't make separate boxes for now */
837:   mesh->startX    = oldMesh->startX;
838:   mesh->endX      = oldMesh->endX;
839:   mesh->startY    = oldMesh->startY;
840:   mesh->endY      = oldMesh->endY;
841:   mesh->locStartX = mesh->startX;
842:   mesh->locEndX   = mesh->endX;
843:   mesh->locStartY = mesh->startY;
844:   mesh->locEndY   = mesh->endY;
845:   mesh->sizeX     = mesh->endX    - mesh->startX;
846:   mesh->sizeY     = mesh->endY    - mesh->startY;
847:   mesh->locSizeX  = mesh->locEndX - mesh->locStartX;
848:   mesh->locSizeY  = mesh->locEndY - mesh->locStartY;

850:   /* Generate a refined mesh */
851:   if (rank == 0) {
852:     /* Generate the refined mesh and check whether we should create midnodes */
853:     if (oldMesh->numCorners == 3) {
854:       PetscStrcpy(cmd_line, "qcenzQ");
855:     } else if (oldMesh->numCorners == 6) {
856:       PetscStrcpy(cmd_line, "qcenzo2Q");
857:     } else {
858:       SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", oldMesh->numCorners);
859:     }
860:     PetscOptionsGetString("ref_", "-triangle_cmd", cmd_line, 251, &opt);
861:     PetscStrcat(cmd_line, "pra");
862:     triangulate(cmd_line, &in, &out, PETSC_NULL);

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

867:     /* Get rid of extra information */
868:     PetscFree(out.segmentlist);
869:     PetscFree(out.segmentmarkerlist);
870:   }
871:   MeshFinalizeRefineInput_Triangle(oldMesh, &in);

873:   /* Communicate values */
874:   MeshDistribute_Triangle(oldMesh, &out);

876:   /* Store the information */
877:   mesh->numBd       = oldMesh->numBd;
878:   mesh->numNodes    = out.numberofpoints;
879:   mesh->numEdges    = out.numberofedges;
880:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
881:   mesh->numFaces    = out.numberoftriangles;
882:   mesh->numCorners  = out.numberofcorners;
883:   tri->nodes       = out.pointlist;
884:   tri->markers     = out.pointmarkerlist;
885:   tri->edges       = out.edgelist;
886:   tri->edgemarkers = out.edgemarkerlist;
887:   tri->faces       = out.trianglelist;
888:   tri->neighbors   = out.neighborlist;
889:   PetscLogObjectMemory(mesh, (mesh->numNodes*2 + mesh->numNodes) * sizeof(double));
890:   PetscLogObjectMemory(mesh, (mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners + mesh->numFaces*3) * sizeof(int));
891:   tri->nodesOld    = PETSC_NULL;

893:   /* Copy holes from previous mesh */
894:   mesh->numHoles    = oldMesh->numHoles;
895:   mesh->holes       = PETSC_NULL;
896:   if (mesh->numHoles > 0) {
897:     PetscMalloc(mesh->numHoles*2 * sizeof(double), &mesh->holes);
898:     PetscMemcpy(mesh->holes, oldMesh->holes, mesh->numHoles*2 * sizeof(double));
899:     PetscLogObjectMemory(mesh, oldMesh->numHoles*2 * sizeof(double));
900:   }

902:   /* Calculate maximum degree of vertices */
903:   MeshSetupSupport_Triangular_2D(mesh);

905:   /* Construct derived and boundary information */
906:   MeshSetupBoundary_Triangular_2D(mesh);

908: #ifdef PETSC_USE_BOPT_g
909:   /* Check mesh integrity */
910:   MeshDebug_Triangular_2D(mesh, PETSC_FALSE);
911: #endif

913:   /* Reorder nodes before distributing mesh */
914:   PetscOptionsHasName(oldMesh->prefix, "-mesh_reorder", &opt);
915:   if (opt == PETSC_TRUE) {
916:     /* MUST FIX: Since we have duplicated the whole mesh, we may impersonate a serial mesh */
917:     comm       = mesh->comm;
918:     mesh->comm = PETSC_COMM_SELF;
919:     ierr       = MeshGetOrdering(mesh, MESH_ORDER_TRIANGULAR_2D_RCM, &mesh->nodeOrdering);
920:     mesh->comm = comm;

922:     /* Permute arrays implicitly numbered by node numbers */
923:     AOApplicationToPetscPermuteReal(mesh->nodeOrdering, 2, tri->nodes);
924:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->markers);
925:     AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->degrees);
926:     /* Renumber arrays dependent on the canonical node numbering */
927:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numEdges*2,                   tri->edges);
928:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numFaces*oldMesh->numCorners, tri->faces);
929:     AOApplicationToPetsc(mesh->nodeOrdering, mesh->numBdNodes,                   tri->bdNodes);
930:   }

932:   /* Initialize partition */
933:   MeshPartition(mesh);

935:   /* Reset the midnodes */
936:   MeshResetNodes(mesh, PETSC_TRUE);
937:   return(0);
938: }

940: #endif /* PETSC_HAVE_TRIANGLE */