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

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

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

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

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

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

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

109: #if 0
110: #ifdef PETSC_USE_BOPT_g
111:   {
112:     int neighbor, neighbor2, row2;

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

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

149:   return(0);
150: }

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

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

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

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

219:           /* Account for constrained nodes: link to anchor */
220:           if ((constrain == PETSC_TRUE) && (markers[sNode] < 0)) {
221:             ierr  = MeshGetBoundaryIndex(mesh, markers[sNode], &sBd);
222:             sNode = bdNodes[bdBegin[sBd]];
223:           }

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

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

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

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

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

316:           /* Account for constrained nodes: link to anchor */
317:           if ((constrain == PETSC_TRUE) && (markers[sNode] < 0)) {
318:             ierr  = MeshGetBoundaryIndex(mesh, markers[sNode], &sBd);
319:             sNode = bdNodes[bdBegin[sBd]];
320:           }

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

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

382:   *numV          = numNodes;
383:   *numE          = nodeOffsets[numNodes]/2;
384:   *vertOffsets   = nodeOffsets;
385:   *vertNeighbors = nodeNeighbors;
386:   return(0);
387: }

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

412:   /* Get the size of the neighbors array */
413:   for(elem = p->firstElement[p->rank], size = 0; elem < p->firstElement[p->rank+1]; elem++) {
414:     for(neighbor = 0; neighbor < numNeighbors; neighbor++) {
415:       if (neighbors[elem*numNeighbors + neighbor] >= 0) size++;
416:     }
417:   }

419:   /* Allocate memory */
420:   numLocElements = p->firstElement[p->rank+1] - p->firstElement[p->rank];
421:   PetscMalloc((numLocElements+1) * sizeof(int), &eO);
422:   PetscMalloc(size               * sizeof(int), &eN);

424:   if ((weight != 0) && (edgeWeights != PETSC_NULL)) {
425:     PetscMalloc(size * sizeof(int), &eW);
426:   } else {
427:     eW   = PETSC_NULL;
428:   }

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

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

482:   *elemOffsets   = eO;
483:   *elemNeighbors = eN;
484:   *edgeWeights   = eW;

486:   return(0);
487: }

489: /*-------------------------------------- Replacement Backend Functions ----------------------------------------------*/
490: #undef  __FUNCT__
492: int PartitionCreate_CSR(Partition p)
493: {
494:   Partition_Triangular_2D *q        = (Partition_Triangular_2D *) p->data;
495:   Mesh                     mesh     = p->mesh;
496:   int                      numProcs = p->numProcs;
497:   int                      proc;
498:   int                      ierr;

501:   p->numLocElements     = mesh->numFaces;
502:   p->numOverlapElements = p->numLocElements;
503:   q->numLocNodes        = mesh->numNodes;
504:   q->numOverlapNodes    = q->numLocNodes;
505:   q->numLocEdges        = mesh->numEdges;
506:   q->numLocBdNodes      = mesh->numBdNodes;
507:   q->numOverlapBdNodes  = q->numLocBdNodes;

509:   MPI_Allgather(&p->numLocElements, 1, MPI_INT, &p->firstElement[1], 1, MPI_INT, p->comm);
510:   MPI_Allgather(&q->numLocNodes,    1, MPI_INT, &q->firstNode[1],    1, MPI_INT, p->comm);
511:   MPI_Allgather(&q->numLocEdges,    1, MPI_INT, &q->firstEdge[1],    1, MPI_INT, p->comm);
512:   MPI_Allgather(&q->numLocBdNodes,  1, MPI_INT, &q->firstBdNode[1],  1, MPI_INT, p->comm);
513:   for(proc = 1; proc < numProcs; proc++) {
514:     p->firstElement[proc+1] += p->firstElement[proc];
515:     q->firstNode[proc+1]    += q->firstNode[proc];
516:     q->firstEdge[proc+1]    += q->firstEdge[proc];
517:     q->firstBdNode[proc+1]  += q->firstBdNode[proc];
518:   }

520:   p->numElements = p->firstElement[numProcs];
521:   q->numNodes    = q->firstNode[numProcs];
522:   q->numEdges    = q->firstEdge[numProcs];
523:   q->numBdNodes  = q->firstBdNode[numProcs];
524:   return(0);
525: }

527: #undef  __FUNCT__
529: int MeshPartition_Triangular_2D_CSR(Mesh mesh)
530: {

534:   PartitionCreate(mesh, &mesh->part);
535:   PetscObjectComposeFunction((PetscObject) mesh->part, "PartitionTriangular2D_Create_C", "PartitionCreate_CSR",
536:                                     (void (*)(void)) PartitionCreate_CSR);
537: 
538:   PartitionSetType(mesh->part, PARTITION_TRIANGULAR_2D);
539:   mesh->partitioned = 1;
540:   PetscFunctionReturn(ierr);
541: }

543: /*-------------------------------------------- Standard Interface ---------------------------------------------------*/
544: #undef  __FUNCT__
546: /* MeshCreate_CSR
547:   This function creates a 2D unstructured mesh using an existing CSR representation.

549:   Collective on Mesh

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

554:   Input Parameters from bdCtx:
555: + numBD       - The number of closed boundaries
556: . numVertices - The number of mesh nodes
557: . vertices    - The (x,y) coordinates of the mesh nodes
558: . markers     - The offset of the adjacency list for each node
559: . segMarkers  - The adjacency list of the mesh
560: . numSegments - The number of elements in the mesh
561: . segments    - The nodes in each element

563:   Output Parameter:
564: . mesh      - The new mesh

566:   Level: developer

568: .keywords mesh, CSR
569: .seealso MeshCreate_Triangle
570: */
571: int MeshCreate_CSR(MeshBoundary2D *bdCtx, int numCorners, Mesh mesh)
572: {
573:   Mesh_Triangular *tri      = (Mesh_Triangular *) mesh->data;
574:   int              numNodes = bdCtx->numVertices;
575:   double          *nodes    = bdCtx->vertices;
576:   int             *offsets  = bdCtx->markers;
577:   int             *adj      = bdCtx->segMarkers;
578:   int              numEdges = offsets[numNodes];
579:   int              numFaces = bdCtx->numSegments;
580:   int             *faces    = bdCtx->segments;
581:   int              size, rank;
582:   int              node, adjNode, edge;
583:   int              ierr;


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

591:   /* Setup function table */
592:   mesh->ops->partition = MeshPartition_Triangular_2D_CSR;

594:   /* Allocate arrays */
595:   tri->nodes         = PETSC_NULL;
596:   tri->markers       = PETSC_NULL;
597:   tri->edges         = PETSC_NULL;
598:   tri->edgemarkers   = PETSC_NULL;
599:   tri->faces         = PETSC_NULL;
600:   if (numNodes > 0) {
601:     PetscMalloc(numNodes*2 * sizeof(double), &tri->nodes);
602:     PetscMalloc(numNodes   * sizeof(int),    &tri->markers);
603:   }
604:   if (numEdges > 0) {
605:     PetscMalloc(numEdges*2 * sizeof(int),    &tri->edges);
606:     PetscMalloc(numEdges   * sizeof(int),    &tri->edgemarkers);
607:   }
608:   if (numFaces > 0) {
609:     PetscMalloc(numFaces*3 * sizeof(int),    &tri->faces);
610:   }

612:   /* Remove redundant and self edges */
613:   for(node = 0, edge = 0; node < numNodes; node++) {
614:     for(adjNode = offsets[node]; adjNode < offsets[node+1]; adjNode++) {
615:       if (adj[adjNode] > node) {
616:         /* Could sort here */
617:         tri->edges[edge*2]   = node;
618:         tri->edges[edge*2+1] = adj[adjNode];
619:         edge++;
620:       }
621:     }
622:   }
623:   numEdges = edge;

625:   /* Store the information */
626:   mesh->numBd       = bdCtx->numBd;
627:   mesh->numNodes    = numNodes;
628:   mesh->numEdges    = numEdges;
629:   mesh->numVertices = (mesh->numNodes <= mesh->numEdges ? mesh->numNodes : mesh->numNodes - mesh->numEdges);
630:   mesh->numFaces    = numFaces;
631:   mesh->numCorners  = 3;
632:   mesh->numHoles    = 0;
633:   PetscMemcpy(tri->nodes, nodes, numNodes*2 * sizeof(double));
634:   PetscMemcpy(tri->faces, faces, numFaces*3 * sizeof(int));
635:   PetscMemzero(tri->markers,    numNodes   * sizeof(int));
636:   PetscMemzero(tri->edgemarkers, numEdges * sizeof(int));
637:   tri->neighbors   = PETSC_NULL;
638:   mesh->holes      = PETSC_NULL;
639:   PetscLogObjectMemory(mesh, (mesh->numNodes*2) * sizeof(double));
640:   PetscLogObjectMemory(mesh, (mesh->numNodes + mesh->numEdges*2 + mesh->numEdges + mesh->numFaces*mesh->numCorners) * sizeof(int));

642:   return(0);
643: }

645: #undef  __FUNCT__
647: /*
648:   MeshRefine_CSR - This function refines a two dimensional unstructured mesh
649:   with a CSR representation using area constraints.
650:   
651:   Collective on Mesh
652:   
653:   Input Parameters:
654: + oldMesh - The mesh begin refined
655: - area    - A function which gives an area constraint when evaluated inside an element
656:   
657:   Output Parameter:
658: . mesh    - The refined mesh
659:   
660:   Level: developer 
661:   
662: .keywords: mesh, CSR, refine 
663: .seealso: MeshRefine(), MeshCoarsen(), MeshReform()
664: */
665: int MeshRefine_CSR(Mesh oldMesh, PointFunction area, Mesh mesh)
666: {
667:   SETERRQ(PETSC_ERR_SUP, "Not yet implemented");
668: }