Actual source code: tri2dCSR.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: tri2dCSR.c,v 1.5 2001/01/10 23:08:38 knepley Exp $";
  3: #endif

  5: /* CSR format for 2d triangular grids */
  6: #include "src/mesh/impls/triangular/2d/2dimpl.h" /*I "mesh.h" I*/
  7: #include "tri2dCSR.h"

  9: int MeshCreateLocalCSR_Triangular_2D(Mesh mesh, int *numV, int *numE, int **vertOffsets, int **vertNeighbors,
 10:                                      double **vertCoords, int numBC, int *bcNodes, PetscTruth symmetric)
 11: {
 12:   Mesh_Triangular *tri         = (Mesh_Triangular *) mesh->data;
 13:   int              numEdges    = mesh->numEdges;
 14:   int             *edges       = tri->edges;
 15:   int              numVertices = mesh->numVertices;
 16:   int              numNodes    = mesh->numNodes;
 17:   double          *nodes       = tri->nodes;
 18:   int             *degrees     = tri->degrees;
 19:   int             *vO;
 20:   int             *vN;
 21:   double          *vC;
 22:   int              size, count, locCount, bcSeen, nBCSeen;
 23: #ifdef PETSC_USE_BOPT_g
 24:   int              globalRows;
 25: #endif
 26:   int              node, nNode, edge, row, bc;
 27:   int              ierr;

 30:   /* Get the size of the neighbors array -- Remember that midnodes have degree 0 */
 31:   for(node = 0, size = 0; node < numNodes; node++) {
 32:     size += degrees[node]+1;
 33:   }

 35:   /* Allocate memory */
 36:   PetscMalloc((numVertices-numBC+1) * sizeof(int),    &vO);
 37:   PetscMalloc(size                  * sizeof(int),    &vN);
 38:   PetscMalloc(numVertices*2         * sizeof(double), &vC);

 40:   /* Fill the neighbors array */
 41:   vO[0] = 0;
 42:   for(node = 0, row = 1, count = 0, bcSeen = 0; node < numNodes; node++) {
 43:     /* Only process vertices */
 44:     if (degrees[node] == 0)
 45:       continue;
 46:     /* Ignore constrained nodes */
 47:     for(bc = 0; bc < numBC; bc++)
 48:       if (node == bcNodes[bc])
 49:         break;
 50:     if (bc < numBC) {
 51:       bcSeen++;
 52:       continue;
 53:     }

 55:     /* Include a self edge */
 56:     vN[count++] = node - bcSeen;
 57: #ifdef PETSC_USE_BOPT_g
 58:     if (node - bcSeen < 0) SETERRQ(PETSC_ERR_PLIB, "Invalid node number");
 59: #endif

 61:     for(edge = 0, locCount = 0; (edge < numEdges) && (locCount < degrees[node]); edge++)
 62:     {
 63:       if (edges[edge*2] == node)
 64:         nNode = edges[edge*2+1];
 65:       else if (edges[edge*2+1] == node)
 66:         nNode = edges[edge*2];
 67:       else
 68:         continue;

 70:       locCount++;
 71:       /* Ignore edges which go off processor and check for symmetrized matrix */
 72:       if ((nNode >= numNodes) || ((symmetric == PETSC_TRUE) && (nNode > node)))
 73:         continue;
 74:       /* Ignore constrained nodes and renumber */
 75:       for(bc = 0, nBCSeen = 0; bc < numBC; bc++) {
 76:         if (bcNodes[bc] < 0)
 77:           continue;
 78:         if (nNode == bcNodes[bc])
 79:           break;
 80:         else if (nNode > bcNodes[bc])
 81:           nBCSeen++;
 82:       }
 83:       if (bc < numBC)
 84:         continue;
 85:       vN[count++] = nNode - nBCSeen;
 86: #ifdef PETSC_USE_BOPT_g
 87:       if (nNode - nBCSeen < 0) SETERRQ(PETSC_ERR_PLIB, "Invalid node number");
 88: #endif
 89:     }
 90:     /* Get coordinates */
 91:     vC[(row-1)*2]   = nodes[node*2];
 92:     vC[(row-1)*2+1] = nodes[node*2+1];
 93:     vO[row++]       = count;
 94:   }
 95: #ifdef PETSC_USE_BOPT_g
 96:   MPI_Allreduce(&row, &globalRows, 1, MPI_INT, MPI_SUM, mesh->comm);
 97:   if (globalRows != numVertices - numBC + mesh->part->numProcs) {
 98:     SETERRQ(PETSC_ERR_PLIB, "Invalid number of vertices or bad boundary conditions");
 99:   }
100: #endif
101:   if (count > size) SETERRQ(PETSC_ERR_PLIB, "Invalid connectivity");
102:   PetscLogInfo(mesh, "Local CSR Graph: %d vertices, %d edgesn", row-1, count-(row-1));
103:   *numV = row-1;
104:   *numE = count-(row-1);

106: #if 0
107: #ifdef PETSC_USE_BOPT_g
108:   {
109:     int neighbor, neighbor2, row2;

111:   /* Check local symmetry of the graph */
112:   for(row = 0; row < numNodes; row++)
113:     for(neighbor = vO[row]; neighbor < vO[row+1]; neighbor++)
114:     {
115:       row2 = vN[neighbor];
116:       if ((row2 < 0) || (row2 >= numNodes)) {
117:         SETERRQ(PETSC_ERR_PLIB, "Invalid local graph");
118:       } else {
119:         /* Check for companion edge */
120:         for(neighbor2 = vO[row2]; neighbor2 < vO[row2+1]; neighbor2++)
121:           if (vN[neighbor2] == row)
122:             break;
123:         if (neighbor2 == vO[row2+1]) SETERRQ(PETSC_ERR_PLIB, "Nonsymmetric local graph");
124:       }
125:     }
126:   }
127: #endif
128: #endif

130:   if (vertOffsets   != PETSC_NULL) {
131:     *vertOffsets   = vO;
132:   } else {
133:     PetscFree(vO);
134:   }
135:   if (vertNeighbors != PETSC_NULL) {
136:     *vertNeighbors = vN;
137:   } else {
138:     PetscFree(vN);
139:   }
140:   if (vertCoords    != PETSC_NULL) {
141:     *vertCoords    = vC;
142:   } else {
143:     PetscFree(vC);
144:   }

146:   return(0);
147: }

149: int MeshCreateFullCSR_Triangular_2D(Mesh mesh, PetscTruth constrain, int *numV, int *numE, int **vertOffsets, int **vertNeighbors)
150: {
151:   Mesh_Triangular *tri         = (Mesh_Triangular *) mesh->data;
152:   int              numBd       = mesh->numBd;
153:   int              numElements = mesh->numFaces;
154:   int             *elements    = tri->faces;
155:   int              numCorners  = mesh->numCorners;
156:   int              numNodes    = mesh->numNodes;
157:   int             *markers     = tri->markers;
158:   int             *bdMarkers   = tri->bdMarkers;
159:   int             *bdBegin     = tri->bdBegin;
160:   int             *bdNodes     = tri->bdNodes;
161:   int             *nodeOffsets, *nodeNeighbors;
162:   int            **constNeighbors;
163:   int             *nodeDone, *neighbors, *constCount;
164:   int             *support;
165:   int              degree;
166:   PetscTruth       duplicate;
167:   int              maxNeighbors, size;
168:   int              bd, elem, corner, neighbor, node, anchorNode, sElem, sCorner, sNode, sBd, nNode, bdNode, count = 0;
169:   int              ierr;

172:   /* Allocate adjacency structures */
173:   maxNeighbors  = mesh->maxDegree*numCorners;
174:   PetscMalloc(numNodes     * sizeof(int), &nodeDone);
175:   PetscMalloc(maxNeighbors * sizeof(int), &nodeNeighbors);
176:   PetscMalloc((numNodes+1) * sizeof(int), &nodeOffsets);
177:   PetscMemzero(nodeOffsets, (numNodes+1) * sizeof(int));
178:   PetscMemzero(nodeDone,     numNodes    * sizeof(int));
179:   /* Allocate constraint structures */
180:   if (constrain == PETSC_TRUE) {
181:     PetscMalloc(numBd * sizeof(int *), &constNeighbors);
182:     for(bd = 0; bd < numBd; bd++) {
183:       if (bdMarkers[bd] < 0) {
184:         PetscMalloc(maxNeighbors*(bdBegin[bd+1] - bdBegin[bd]) * sizeof(int), &constNeighbors[bd]);
185:       }
186:     }
187:   }
188:   /* Calculate the number of neighbors for each node */
189:   for(elem = 0; elem < numElements; elem++) {
190:     for(corner = 0; corner < numCorners; corner++) {
191:       node = elements[elem*numCorners+corner];
192:       if (nodeDone[node]) continue;
193:       nodeDone[node] = 1;

195:       /* Get array holding list of neighboring nodes */
196:       if ((constrain == PETSC_TRUE) && (markers[node] < 0)) {
197:         /* Constrained node: We maintain a list of neighbors seen */
198:         MeshGetBoundaryIndex(mesh, markers[node], &bd);
199:         neighbors  = constNeighbors[bd];
200:         /* We let the first constrained node be a representative for all of them */
201:         anchorNode = bdNodes[bdBegin[bd]];
202:       } else {
203:         /* Normal node: Just use temporary storage since we only look at it once */
204:         neighbors  = nodeNeighbors;
205:         anchorNode = node;
206:       }

208:       /* Loop over nodes on each element in the support of the node */
209:       MeshGetNodeSupport(mesh, node, elem, &degree, &support);
210:       for(sElem = 0, count = nodeOffsets[anchorNode+1]; sElem < degree; sElem++) {
211:         for(sCorner = 0; sCorner < numCorners; sCorner++) {
212:           sNode = elements[support[sElem]*numCorners+sCorner];

214:           /* Account for constrained nodes: link to anchor */
215:           if ((constrain == PETSC_TRUE) && (markers[sNode] < 0)) {
216:             ierr  = MeshGetBoundaryIndex(mesh, markers[sNode], &sBd);
217:             sNode = bdNodes[bdBegin[sBd]];
218:           }

220:           /* Check for duplicate node and insert in order */
221:           for(neighbor = count-1, duplicate = PETSC_FALSE; neighbor >= 0; neighbor--) {
222:             nNode = neighbors[neighbor];
223:             if (sNode > nNode) {
224:               break;
225:             } else if (sNode == nNode) {
226:               duplicate = PETSC_TRUE;
227:               break;
228:             }
229:             neighbors[neighbor+1] = nNode;
230:           }
231:           if (duplicate == PETSC_FALSE) {
232:             neighbors[neighbor+1] = sNode;
233:             count++;
234:           } else {
235:             PetscMemmove(&neighbors[neighbor+1], &neighbors[neighbor+2], (count-1-neighbor) * sizeof(int));
236: 
237:           }
238: #ifdef PETSC_USE_BOPT_g
239:           if ((constrain == PETSC_TRUE) && (markers[node] < 0)) {
240:             if (count >= maxNeighbors*(bdBegin[bd+1] - bdBegin[bd]))
241:               SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many neighboring nodes %d", count);
242:           } else {
243:             if (count >= maxNeighbors)
244:               SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many neighboring nodes %d", count);
245:           }
246: #endif
247:         }
248:       }
249:       MeshRestoreNodeSupport(mesh, node, elem, &degree, &support);

251:       /* Record number of neighbors */
252:       nodeOffsets[anchorNode+1] = count;
253:     }
254:   }
255:   PetscFree(nodeNeighbors);
256:   /* Connect all the nodes on a given boundary */
257:   if (constrain == PETSC_TRUE) {
258:     for(bd = 0; bd < numBd; bd++) {
259:       if (bdMarkers[bd] < 0) {
260:         anchorNode = bdNodes[bdBegin[bd]];
261:         for(bdNode = bdBegin[bd]+1; bdNode < bdBegin[bd+1]; bdNode++)
262:           nodeOffsets[bdNodes[bdNode]+1] = nodeOffsets[anchorNode+1];
263:       }
264:     }
265:   }
266:   /* Do prefix sums */
267:   nodeOffsets[0] = 0;
268:   for(node = 1; node <= numNodes; node++) {
269:     nodeOffsets[node] += nodeOffsets[node-1];
270:   }
271:   /* Cleanup */
272:   if (constrain == PETSC_TRUE) {
273:     for(bd = 0; bd < numBd; bd++) {
274:       if (bdMarkers[bd] < 0) {
275:         PetscFree(constNeighbors[bd]);
276:       }
277:     }
278:     PetscFree(constNeighbors);
279:   }

281:   /* Calculate adjacency list */
282:   PetscMalloc(nodeOffsets[numNodes] * sizeof(int), &nodeNeighbors);
283:   PetscMalloc(numBd                 * sizeof(int), &constCount);
284:   PetscMemzero(nodeDone,   numNodes * sizeof(int));
285:   PetscMemzero(constCount, numBd    * sizeof(int));
286:   for(elem = 0; elem < numElements; elem++) {
287:     for(corner = 0; corner < numCorners; corner++) {
288:       node = elements[elem*numCorners+corner];
289:       if (nodeDone[node]) continue;
290:       nodeDone[node] = 1;

292:       /* Get array holding list of neighboring nodes */
293:       if ((constrain == PETSC_TRUE) && (markers[node] < 0)) {
294:         /* We let the first constrained node be a representative for all of them */
295:         MeshGetBoundaryIndex(mesh, markers[node], &bd);
296:         anchorNode = bdNodes[bdBegin[bd]];
297:         count      = constCount[bd];
298:       } else {
299:         anchorNode = node;
300:         count      = 0;
301:       }
302:       neighbors    = &nodeNeighbors[nodeOffsets[anchorNode]];

304:       /* Loop over nodes on each element in the support of the node */
305:       MeshGetNodeSupport(mesh, node, elem, &degree, &support);
306:       for(sElem = 0; sElem < degree; sElem++) {
307:         for(sCorner = 0; sCorner < numCorners; sCorner++) {
308:           if (count >= nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode]) break;
309:           sNode = elements[support[sElem]*numCorners+sCorner];

311:           /* Account for constrained nodes: link to anchor */
312:           if ((constrain == PETSC_TRUE) && (markers[sNode] < 0)) {
313:             ierr  = MeshGetBoundaryIndex(mesh, markers[sNode], &sBd);
314:             sNode = bdNodes[bdBegin[sBd]];
315:           }

317:           /* Check for duplicate node and insert in order */
318:           for(neighbor = count-1, duplicate = PETSC_FALSE; neighbor >= 0; neighbor--) {
319:             nNode = neighbors[neighbor];
320:             if (sNode > nNode) {
321:               break;
322:             } else if (sNode == nNode) {
323:               duplicate = PETSC_TRUE;
324:               break;
325:             }
326:             neighbors[neighbor+1] = nNode;
327:           }
328:           if (duplicate == PETSC_FALSE) {
329:             neighbors[neighbor+1] = sNode;
330:             count++;
331:           } else {
332:             PetscMemmove(&neighbors[neighbor+1], &neighbors[neighbor+2], (count-1-neighbor) * sizeof(int));
333: 
334:           }
335:         }
336:       }
337:       MeshRestoreNodeSupport(mesh, node, elem, &degree, &support);

339:       if ((constrain == PETSC_FALSE) || (markers[node] >= 0)) {
340:         if (count != nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode]) {
341:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of adjacent nodes %d should be %d",
342:                    count, nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode]);
343:         }
344:       } else {
345:         constCount[bd] = count;
346:       }
347:     }
348:   }
349:   /* Handle constrained boundaries */
350:   if (constrain == PETSC_TRUE) {
351:     for(bd = 0; bd < numBd; bd++) {
352:       if (bdMarkers[bd] < 0) {
353:         anchorNode = bdNodes[bdBegin[bd]];
354:         /* Check adjacency list */
355:         if (constCount[bd] != nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode]) {
356:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of adjacent nodes %d should be %d",
357:                    count, nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode]);
358:         }
359:         /* Connect all the nodes on a given boundary */
360:         for(bdNode = bdBegin[bd]+1; bdNode < bdBegin[bd+1]; bdNode++) {
361:           node = bdNodes[bdNode];
362:           size = nodeOffsets[anchorNode+1] - nodeOffsets[anchorNode];
363:           if (nodeOffsets[node+1] - nodeOffsets[node] != size) {
364:             SETERRQ2(PETSC_ERR_PLIB, "Invalid number of adjacent nodes %d should be %d",
365:                      nodeOffsets[node+1] - nodeOffsets[node], size);
366:           }
367:           PetscMemcpy(&nodeNeighbors[nodeOffsets[node]], &nodeNeighbors[nodeOffsets[anchorNode]], size*sizeof(int));
368: 
369:         }
370:       }
371:     }
372:   }
373:   /* Cleanup */
374:   PetscFree(nodeDone);
375:   PetscFree(constCount);

377:   *numV          = numNodes;
378:   *numE          = nodeOffsets[numNodes]/2;
379:   *vertOffsets   = nodeOffsets;
380:   *vertNeighbors = nodeNeighbors;
381:   return(0);
382: }

384: int MeshCreateDualCSR_Triangular_2D(Mesh mesh, int **elemOffsets, int **elemNeighbors, int **edgeWeights, int weight)
385: {
386:   Mesh_Triangular *tri          = (Mesh_Triangular *) mesh->data;
387:   Partition        p            = mesh->part;
388:   int              numNeighbors = 3;
389:   int              numCorners   = mesh->numCorners;
390:   int             *neighbors    = tri->neighbors;
391:   int             *elements     = tri->faces;
392:   int             *markers      = tri->markers;
393:   int             *eO;
394:   int             *eN;
395:   int             *eW;
396:   int              numLocElements;
397:   int              size, count;
398:   int              elem, neighbor, corner, row;
399: #ifdef PETSC_USE_BOPT_g
400:   int              neighbor2, row2;
401: #endif
402:   int              ierr;

405:   /* Get the size of the neighbors array */
406:   for(elem = p->firstElement[p->rank], size = 0; elem < p->firstElement[p->rank+1]; elem++) {
407:     for(neighbor = 0; neighbor < numNeighbors; neighbor++) {
408:       if (neighbors[elem*numNeighbors + neighbor] >= 0) size++;
409:     }
410:   }

412:   /* Allocate memory */
413:   numLocElements = p->firstElement[p->rank+1] - p->firstElement[p->rank];
414:   PetscMalloc((numLocElements+1) * sizeof(int), &eO);
415:   PetscMalloc(size               * sizeof(int), &eN);

417:   if ((weight != 0) && (edgeWeights != PETSC_NULL)) {
418:     PetscMalloc(size * sizeof(int), &eW);
419:   } else {
420:     eW   = PETSC_NULL;
421:   }

423:   /* Fill the neighbors array */
424:   eO[0] = 0;
425:   for(elem = p->firstElement[p->rank], row = 1, count = 0; elem < p->firstElement[p->rank+1]; elem++, row++) {
426:     for(neighbor = 0; neighbor < numNeighbors; neighbor++) {
427:       if (neighbors[elem*numNeighbors + neighbor] >= 0) {
428:         eN[count] = neighbors[elem*numNeighbors + neighbor];
429:         if (eW != PETSC_NULL) {
430:           /* Check element for a node on an inner boundary */
431:           for(corner = 0; corner < 3; corner++) {
432:             if (markers[elements[elem*numCorners+corner]] < 0) break;
433:           }
434:           if (corner < 3) {
435:             /* Check neighbor for a node on an inner boundary */
436:             for(corner = 0; corner < 3; corner++) {
437:               if (markers[elements[neighbors[elem*numNeighbors + neighbor]*numCorners+corner]] < 0) break;
438:             }
439:             if (corner < 3) {
440:               eW[count] = weight;
441:             } else {
442:               eW[count] = 0;
443:             }
444:           } else {
445:             eW[count] = 0;
446:           }
447:         }
448:         count++;
449:       }
450:     }
451:     eO[row] = count;
452:   }
453:   if (count != size) SETERRQ(PETSC_ERR_PLIB, "Invalid adjacency matrix");
454:   PetscLogInfo(mesh, "Dual: %d elements, %d edgesn", numLocElements, size);

456: #ifdef PETSC_USE_BOPT_g
457:   /* Check local symmetry of the dual graph */
458:   for(row = 0; row < numLocElements; row++) {
459:     for(neighbor = eO[row]; neighbor < eO[row+1]; neighbor++) {
460:       row2 = eN[neighbor] - p->firstElement[p->rank];
461:       if ((row2 < 0) || (row2 >= numLocElements)) {
462:         /* PetscSynchronizedPrintf(p->comm, "[%3d]Cut: %6d --> %6dn", p->rank, row + p->firstElement[p->rank], eN[neighbor]); */
463:       } else {
464:         /* Check for companion edge */
465:         for(neighbor2 = eO[row2]; neighbor2 < eO[row2+1]; neighbor2++) {
466:           if (eN[neighbor2] - p->firstElement[p->rank] == row) break;
467:         }
468:         if (neighbor2 == eO[row2+1]) SETERRQ(PETSC_ERR_PLIB, "Nonsymmetric dual graph");
469:       }
470:     }
471:   }
472:   /* PetscSynchronizedFlush(p->comm); */
473: #endif

475:   *elemOffsets   = eO;
476:   *elemNeighbors = eN;
477:   *edgeWeights   = eW;

479:   return(0);
480: }

482: /*-------------------------------------- Replacement Backend Functions ----------------------------------------------*/
483: int PartitionCreate_CSR(Partition p)
484: {
485:   Partition_Triangular_2D *q        = (Partition_Triangular_2D *) p->data;
486:   Mesh                     mesh     = p->mesh;
487:   int                      numProcs = p->numProcs;
488:   int                      proc;
489:   int                      ierr;

492:   p->numLocElements     = mesh->numFaces;
493:   p->numOverlapElements = p->numLocElements;
494:   q->numLocNodes        = mesh->numNodes;
495:   q->numOverlapNodes    = q->numLocNodes;
496:   q->numLocEdges        = mesh->numEdges;
497:   q->numLocBdNodes      = mesh->numBdNodes;
498:   q->numOverlapBdNodes  = q->numLocBdNodes;

500:   MPI_Allgather(&p->numLocElements, 1, MPI_INT, &p->firstElement[1], 1, MPI_INT, p->comm);
501:   MPI_Allgather(&q->numLocNodes,    1, MPI_INT, &q->firstNode[1],    1, MPI_INT, p->comm);
502:   MPI_Allgather(&q->numLocEdges,    1, MPI_INT, &q->firstEdge[1],    1, MPI_INT, p->comm);
503:   MPI_Allgather(&q->numLocBdNodes,  1, MPI_INT, &q->firstBdNode[1],  1, MPI_INT, p->comm);
504:   for(proc = 1; proc < numProcs; proc++) {
505:     p->firstElement[proc+1] += p->firstElement[proc];
506:     q->firstNode[proc+1]    += q->firstNode[proc];
507:     q->firstEdge[proc+1]    += q->firstEdge[proc];
508:     q->firstBdNode[proc+1]  += q->firstBdNode[proc];
509:   }

511:   p->numElements = p->firstElement[numProcs];
512:   q->numNodes    = q->firstNode[numProcs];
513:   q->numEdges    = q->firstEdge[numProcs];
514:   q->numBdNodes  = q->firstBdNode[numProcs];
515:   return(0);
516: }

518: int MeshPartition_Triangular_2D_CSR(Mesh mesh)
519: {

523:   PartitionCreate(mesh, &mesh->part);
524:   PetscObjectComposeFunction((PetscObject) mesh->part, "PartitionTriangular2D_Create_C", "PartitionCreate_CSR",
525:                                     (void (*)(void)) PartitionCreate_CSR);
526: 
527:   PartitionSetType(mesh->part, PARTITION_TRIANGULAR_2D);
528:   mesh->partitioned = 1;
529:   PetscFunctionReturn(ierr);
530: }

532: /*-------------------------------------------- Standard Interface ---------------------------------------------------*/
533: /* MeshCreate_CSR
534:   This function creates a 2D unstructured mesh using an existing CSR representation.

536:   Collective on Mesh

538:   Input Parameter:
539: . numCorners  - The number of nodes in an element, here 3

541:   Input Parameters from bdCtx:
542: + numBD       - The number of closed boundaries
543: . numVertices - The number of mesh nodes
544: . vertices    - The (x,y) coordinates of the mesh nodes
545: . markers     - The offset of the adjacency list for each node
546: . segMarkers  - The adjacency list of the mesh
547: . numSegments - The number of elements in the mesh
548: . segments    - The nodes in each element

550:   Output Parameter:
551: . mesh      - The new mesh

553:   Level: developer

555: .keywords mesh, CSR
556: .seealso MeshCreate_Triangle
557: */
558: int MeshCreate_CSR(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
559: {
560:   Mesh_Triangular *tri      = (Mesh_Triangular *) mesh->data;
561:   int              numNodes = bdCtx->numVertices;
562:   double          *nodes    = bdCtx->vertices;
563:   int             *offsets  = bdCtx->markers;
564:   int             *adj      = bdCtx->segMarkers;
565:   int              numEdges = offsets[numNodes];
566:   int              numFaces = bdCtx->numSegments;
567:   int             *faces    = bdCtx->segments;
568:   int              size, rank;
569:   int              node, adjNode, edge;
570:   int              ierr;


574:   MPI_Comm_size(mesh->comm, &size);
575:   MPI_Comm_rank(mesh->comm, &rank);
576:   if (numCorners != 3) SETERRQ1(PETSC_ERR_SUP, "Each face must have 3 corners, not %d", numCorners);

578:   /* Setup function table */
579:   mesh->ops->partition = MeshPartition_Triangular_2D_CSR;

581:   /* Allocate arrays */
582:   tri->nodes         = PETSC_NULL;
583:   tri->markers       = PETSC_NULL;
584:   tri->edges         = PETSC_NULL;
585:   tri->edgemarkers   = PETSC_NULL;
586:   tri->faces         = PETSC_NULL;
587:   if (numNodes > 0) {
588:     PetscMalloc(numNodes*2 * sizeof(double), &tri->nodes);
589:     PetscMalloc(numNodes   * sizeof(int),    &tri->markers);
590:   }
591:   if (numEdges > 0) {
592:     PetscMalloc(numEdges*2 * sizeof(int),    &tri->edges);
593:     PetscMalloc(numEdges   * sizeof(int),    &tri->edgemarkers);
594:   }
595:   if (numFaces > 0) {
596:     PetscMalloc(numFaces*3 * sizeof(int),    &tri->faces);
597:   }

599:   /* Remove redundant and self edges */
600:   for(node = 0, edge = 0; node < numNodes; node++) {
601:     for(adjNode = offsets[node]; adjNode < offsets[node+1]; adjNode++) {
602:       if (adj[adjNode] > node) {
603:         /* Could sort here */
604:         tri->edges[edge*2]   = node;
605:         tri->edges[edge*2+1] = adj[adjNode];
606:         edge++;
607:       }
608:     }
609:   }
610:   numEdges = edge;

612:   /* Store the information */
613:   mesh->numBd       = bdCtx->numBd;
614:   mesh->numNodes    = numNodes;
615:   mesh->numEdges    = numEdges;
616:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
617:   mesh->numFaces    = numFaces;
618:   mesh->numCorners  = 3;
619:   mesh->numHoles    = 0;
620:   PetscMemcpy(tri->nodes, nodes, numNodes*2 * sizeof(double));
621:   PetscMemcpy(tri->faces, faces, numFaces*3 * sizeof(int));
622:   PetscMemzero(tri->markers,    numNodes   * sizeof(int));
623:   PetscMemzero(tri->edgemarkers, numEdges * sizeof(int));
624:   tri->neighbors   = PETSC_NULL;
625:   mesh->holes      = PETSC_NULL;
626:   PetscLogObjectMemory(mesh, (mesh->numNodes*2) * sizeof(double));
627:   PetscLogObjectMemory(mesh, (mesh->numNodes + mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners) * sizeof(int));

629:   return(0);
630: }

632: /*
633:   MeshRefine_CSR - This function refines a two dimensional unstructured mesh
634:   with a CSR representation using area constraints.
635:   
636:   Collective on Mesh
637:   
638:   Input Parameters:
639: + oldMesh - The mesh begin refined
640: - area    - A function which gives an area constraint when evaluated inside an element
641:   
642:   Output Parameter:
643: . mesh    - The refined mesh
644:   
645:   Level: developer 
646:   
647: .keywords: mesh, CSR, refine 
648: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
649: */
650: int MeshRefine_CSR(Mesh oldMesh, PointFunction area, Mesh mesh)
651: {
652:   SETERRQ(PETSC_ERR_SUP, "Not yet implemented");
653: }