Actual source code: varorder2d.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: varorder2d.c,v 1.8 2000/01/30 18:31:33 huangp Exp $";
  3: #endif

  5: #include "src/grid/gridimpl.h"         /*I "grid.h" I*/
  6: #include "varorder2d.h"

  8: #undef  __FUNCT__
 10: int GridCreateVarOrdering_Triangular_2D(Grid grid, FieldClassMap map, VarOrdering *order) {
 11:   Mesh                  mesh;
 12:   Partition             part;
 13:   PetscConstraintObject constCtx        = grid->constraintCtx;
 14:   int                   numFields       = map->numFields;
 15:   int                  *fields          = map->fields;
 16:   int                   numNodes        = map->numNodes;
 17:   int                   numOverlapNodes = map->numOverlapNodes;
 18:   int                   numGhostNodes   = map->numGhostNodes;
 19:   int                   numClasses      = map->numClasses;
 20:   int                 **fieldClasses    = map->fieldClasses;
 21:   int                  *classes         = map->classes;
 22:   int                  *classSizes      = map->classSizes;
 23:   int                  *localOffsets;
 24:   int                   numNewVars;
 25:   VarOrdering           o;
 26:   /* Ghost variable communication */
 27:   int                  *ghostSendVars;    /* Number of ghost variables on a given processor interior to this domain */
 28:   int                  *sumSendVars;      /* Prefix sums of ghostSendVars */
 29:   int                  *ghostRecvVars;    /* Number of ghost variables on a given processor */
 30:   int                  *sumRecvVars;      /* Prefix sums of ghostRecvVars */
 31:   int                  *displs;           /* Offsets into ghostRecvVars */
 32:   int                   numSendGhostVars; /* The number of ghost variable offsets to send to other processors */
 33:   int                  *sendGhostBuffer;  /* Recv: Global node numbers Send: Offsets of these nodes */
 34:   int                   numProcs, rank;
 35:   int                   proc, f, field, comp, node, locNode, gNode, nclass, var;
 36:   int                   ierr;

 39:   /* Create the ordering */
 40:   PetscHeaderCreate(o, _VarOrdering, int, VAR_ORDER_COOKIE, 0, "VarOrdering", grid->comm, VarOrderingDestroy, 0);
 41:   PetscLogObjectCreate(o);
 42:   PetscObjectCompose((PetscObject) o, "ClassMap", (PetscObject) map);

 44:   /* Allocate memory */
 45:   MPI_Comm_size(grid->comm, &numProcs);
 46:   MPI_Comm_rank(grid->comm, &rank);
 47:   GridGetNumFields(grid, &o->numTotalFields);
 48:   PetscMalloc((numProcs+1)      * sizeof(int),   &o->firstVar);
 49:   PetscMalloc(numOverlapNodes   * sizeof(int),   &o->offsets);
 50:   PetscMalloc(o->numTotalFields * sizeof(int *), &o->localStart);
 51:   PetscLogObjectMemory(o, (numProcs+1 + numOverlapNodes + o->numTotalFields*numClasses)*sizeof(int) + o->numTotalFields*sizeof(int *));
 52:   PetscMemzero(o->localStart, o->numTotalFields * sizeof(int *));
 53:   o->numLocNewVars = 0;
 54:   o->numNewVars    = 0;

 56:   /* Setup domain variable numbering */
 57:   o->offsets[0] = 0;
 58:   for(node = 1; node < numNodes; node++)
 59:     o->offsets[node] = o->offsets[node-1] + classSizes[classes[node-1]];
 60:   o->numLocVars = o->offsets[numNodes-1] + classSizes[classes[numNodes-1]];
 61:   if (map->isConstrained == PETSC_TRUE) {
 62:     (*constCtx->ops->getsize)(constCtx, &o->numLocNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
 63: 
 64:     o->numLocVars += o->numLocNewVars;
 65:   }
 66:   MPI_Allgather(&o->numLocVars, 1, MPI_INT, &o->firstVar[1], 1, MPI_INT, o->comm);
 67:   o->firstVar[0] = 0;
 68:   for(proc = 1; proc <= numProcs; proc++)
 69:     o->firstVar[proc] += o->firstVar[proc-1];
 70:   o->numVars = o->firstVar[numProcs];
 71:   if (map->isConstrained == PETSC_TRUE) {
 72:     (*constCtx->ops->getsize)(constCtx, PETSC_NULL, &o->numNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
 73: 
 74:     MPI_Allreduce(&o->numLocNewVars, &numNewVars, 1, MPI_INT, MPI_SUM, o->comm);
 75:     if (o->numNewVars != numNewVars)
 76:       SETERRQ(PETSC_ERR_PLIB, "Invalid partition of new variables");
 77:   }

 79:   /* Initialize overlap size */
 80:   o->numOverlapVars    = o->numLocVars;
 81:   o->numOverlapNewVars = o->numLocNewVars;

 83:   GridGetMesh(grid, &mesh);
 84:   MeshGetPartition(mesh, &part);
 85:   if (numProcs > 1) {
 86:     /* Map local to global variable numbers */
 87:     for(node = 0; node < numNodes; node++)
 88:       o->offsets[node] += o->firstVar[rank];

 90: #if 0
 91:     GridGhostExchange(o->comm, numGhostNodes, q->ghostNodeProcs, q->ghostNodes, PETSC_INT,
 92:                              q->firstNode, INSERT_VALUES, SCATTER_FORWARD, o->offsets, &o->offsets[numNodes]);
 93: #else
 94:     /* Initialize communication */
 95:     PetscMalloc(numProcs * sizeof(int), &ghostSendVars);
 96:     PetscMalloc(numProcs * sizeof(int), &sumSendVars);
 97:     PetscMalloc(numProcs * sizeof(int), &ghostRecvVars);
 98:     PetscMalloc(numProcs * sizeof(int), &sumRecvVars);
 99:     PetscMalloc(numProcs * sizeof(int), &displs);
100:     PetscMemzero(ghostSendVars, numProcs * sizeof(int));
101:     PetscMemzero(sumSendVars,   numProcs * sizeof(int));
102:     PetscMemzero(ghostRecvVars, numProcs * sizeof(int));
103:     PetscMemzero(sumRecvVars,   numProcs * sizeof(int));
104:     PetscMemzero(displs,        numProcs * sizeof(int));

106:     /* Get number of ghost variables to receive from each processor and size of blocks --
107:          we here assume that classes[] already has ghost node classes in it */
108:     for(node = 0; node < numGhostNodes; node++) {
109:       PartitionGhostToGlobalNodeIndex(part, node, &gNode, &proc);
110:       nclass = classes[numNodes+node];
111:       ghostRecvVars[proc]++;
112:       o->numOverlapVars += classSizes[nclass];
113:     }

115:     /* Get number of constrained ghost variables to receive from each processor and size of blocks */
116:     if (map->isConstrained == PETSC_TRUE) {
117:       (*constCtx->ops->getsize)(constCtx, PETSC_NULL, PETSC_NULL, &o->numOverlapNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
118: 
119:     }
120:     o->numOverlapVars += o->numOverlapNewVars - o->numLocNewVars;

122:     /* Get sizes of ghost variable blocks to send to each processor */
123:     MPI_Alltoall(ghostRecvVars, 1, MPI_INT, ghostSendVars, 1, MPI_INT, o->comm);

125:     /* Calculate offets into the ghost variable receive array */
126:     for(proc = 1; proc < numProcs; proc++) {
127:       sumRecvVars[proc] = sumRecvVars[proc-1] + ghostRecvVars[proc-1];
128:       displs[proc]      = sumRecvVars[proc];
129:     }

131:     /* Calculate offsets into the ghost variable send array */
132:     for(proc = 1; proc < numProcs; proc++)
133:       sumSendVars[proc] = sumSendVars[proc-1] + ghostSendVars[proc-1];

135:     /* Send requests for ghost variable offsets to each processor */
136:     numSendGhostVars = sumSendVars[numProcs-1] + ghostSendVars[numProcs-1];
137:     PetscMalloc(numSendGhostVars * sizeof(int), &sendGhostBuffer);
138:     for(node = 0; node < numGhostNodes; node++) {
139:       PartitionGhostToGlobalNodeIndex(part, node, &gNode, &proc);
140:       o->offsets[numNodes+(displs[proc]++)] = gNode;
141:     }
142:     MPI_Alltoallv(&o->offsets[numNodes],  ghostRecvVars, sumRecvVars, MPI_INT,
143:                          sendGhostBuffer,        ghostSendVars, sumSendVars, MPI_INT, o->comm);
144: 

146:     /* Send ghost variables offsets to each processor */
147:     for(node = 0; node < numSendGhostVars; node++) {
148:       PartitionGlobalToLocalNodeIndex(part, sendGhostBuffer[node], &locNode);
149:       sendGhostBuffer[node] = o->offsets[locNode];
150:     }
151:     MPI_Alltoallv(sendGhostBuffer,       ghostSendVars, sumSendVars, MPI_INT,
152:                          &o->offsets[numNodes], ghostRecvVars, sumRecvVars, MPI_INT, o->comm);
153: 

155:     /* Cleanup */
156:     PetscFree(ghostSendVars);
157:     PetscFree(sumSendVars);
158:     PetscFree(ghostRecvVars);
159:     PetscFree(sumRecvVars);
160:     PetscFree(displs);
161:     PetscFree(sendGhostBuffer);
162: #endif

164:     /* We maintain local offsets for ghost variables, meaning the offsets after the last
165:        interior variable, rather than the offset of the given ghost variable in the global
166:        matrix. */
167:     PetscMalloc(numGhostNodes * sizeof(int), &o->localOffsets);
168:     for(node = 0, var = o->numLocVars; node < numGhostNodes; node++) {
169:       o->localOffsets[node] = var;
170:       nclass = classes[numNodes+node];
171:       var   += classSizes[nclass];
172:     }
173:   }

175:   /* Allocate memory */
176:   PetscMalloc(numClasses * sizeof(int), &localOffsets);
177:   PetscMemzero(localOffsets, numClasses * sizeof(int));

179:   /* Setup local field offsets */
180:   for(f = 0; f < numFields; f++) {
181:     field = fields[f];
182:     PetscMalloc(numClasses * sizeof(int), &o->localStart[field]);
183:     for(nclass = 0; nclass < numClasses; nclass++) {
184:       if (fieldClasses[f][nclass]) {
185:         GridGetFieldComponents(grid, field, &comp);
186:         o->localStart[field][nclass]  = localOffsets[nclass];
187:         localOffsets[nclass]         += comp;
188:       } else {
189:         o->localStart[field][nclass]  = -1;
190:       }
191:     }
192:   }
193:   *order = o;

195:   /* Cleanup */
196:   PetscFree(localOffsets);

198:   return(0);
199: }

201: #undef  __FUNCT__
203: int GridVarOrderingConstrain_Triangular_2D(Grid grid, int constField, PetscConstraintObject constCtx,
204:                                            FieldClassMap constMap, VarOrdering oldOrder, VarOrdering *order)
205: {
206:   Mesh          mesh;
207:   Partition     part;
208:   int           numFields          = constMap->numFields;
209:   int          *fields             = constMap->fields;
210:   int           numNodes           = constMap->numNodes;
211:   int           numOverlapNodes    = constMap->numOverlapNodes;
212:   int           numGhostNodes      = constMap->numGhostNodes;
213:   int           numClasses         = constMap->numClasses;
214:   int          *classes            = constMap->classes;
215:   int          *classMap           = constMap->classMap[constMap->mapSize-1];
216:   int         **localStart         = oldOrder->localStart;
217:   int           numClassesOld;
218:   int           comp;
219:   FieldClassMap map;
220:   VarOrdering   o;
221:   /* Ghost variable communication */
222:   int          *ghostSendVars;    /* Number of ghost variables on a given processor interior to this domain */
223:   int          *sumSendVars;      /* Prefix sums of ghostSendVars */
224:   int          *ghostRecvVars;    /* Number of ghost variables on a given processor */
225:   int          *sumRecvVars;      /* Prefix sums of ghostRecvVars */
226:   int          *displs;           /* Offsets into ghostRecvVars */
227:   int           numSendGhostVars; /* The number of ghost variable offsets to send to other processors */
228:   int          *sendGhostBuffer;  /* Recv: Global node numbers Send: Offsets of these nodes */
229:   int           numProcs, rank;
230:   int           proc, f, field, node, locNode, gNode, nclass, i, var;
231:   int           ierr;

234:   /* Create the ordering */
235:   PetscHeaderCreate(o, _VarOrdering, int, VAR_ORDER_COOKIE, 0, "VarOrdering", grid->comm, VarOrderingDestroy, 0);
236:   PetscLogObjectCreate(o);
237:   PetscObjectCompose((PetscObject) o, "ClassMap", (PetscObject) constMap);

239:   /* Allocate memory */
240:   MPI_Comm_size(grid->comm, &numProcs);
241:   MPI_Comm_rank(grid->comm, &rank);
242:   GridGetFieldComponents(grid, constField, &comp);
243:   o->numTotalFields = oldOrder->numTotalFields;
244:   PetscMalloc((numProcs+1)      * sizeof(int),   &o->firstVar);
245:   PetscMalloc(numOverlapNodes   * sizeof(int),   &o->offsets);
246:   PetscMalloc(o->numTotalFields * sizeof(int *), &o->localStart);
247:   PetscLogObjectMemory(o, (numProcs+1 + numOverlapNodes + o->numTotalFields*numClasses) * sizeof(int) +
248:                        o->numTotalFields * sizeof(int *));
249:   PetscMemzero(o->localStart, o->numTotalFields * sizeof(int *));

251:   /* Setup domain variable numbering */
252:   if (numOverlapNodes < numNodes) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition");
253:   o->offsets[0] = 0;
254:   for(node = 1; node < numNodes; node++)
255:     o->offsets[node]    = o->offsets[node-1] + constMap->classSizes[classes[node-1]];
256:   o->numLocVars = o->offsets[numNodes-1] + constMap->classSizes[classes[numNodes-1]];
257:   (*constCtx->ops->getsize)(constCtx, &o->numLocNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
258: 
259:   o->numLocVars        += o->numLocNewVars;
260:   MPI_Allgather(&o->numLocVars, 1, MPI_INT, &o->firstVar[1], 1, MPI_INT, o->comm);
261:   o->firstVar[0] = 0;
262:   for(proc = 1; proc <= numProcs; proc++)
263:     o->firstVar[proc] += o->firstVar[proc-1];
264:   o->numVars = o->firstVar[numProcs];
265: #ifdef PETSC_USE_BOPT_g
266:   (*constCtx->ops->getsize)(constCtx, PETSC_NULL, &o->numNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
267: 
268:   MPI_Allreduce(&o->numLocNewVars, &i, 1, MPI_INT, MPI_SUM, o->comm);
269:   if (o->numNewVars != i) SETERRQ(PETSC_ERR_PLIB, "Invalid partition of new variables");
270: #endif

272:   /* Initialize variable overlap */
273:   o->numOverlapVars    = o->numLocVars;
274:   o->numOverlapNewVars = o->numLocNewVars;

276:   GridGetMesh(grid, &mesh);
277:   MeshGetPartition(mesh, &part);
278:   if (numProcs > 1) {
279:     /* Map local to global variable numbers */
280:     for(node = 0; node < numNodes; node++)
281:       o->offsets[node] += o->firstVar[rank];

283:     /* Initialize communication */
284:     PetscMalloc(numProcs * sizeof(int), &ghostSendVars);
285:     PetscMalloc(numProcs * sizeof(int), &sumSendVars);
286:     PetscMalloc(numProcs * sizeof(int), &ghostRecvVars);
287:     PetscMalloc(numProcs * sizeof(int), &sumRecvVars);
288:     PetscMalloc(numProcs * sizeof(int), &displs);
289:     PetscMemzero(ghostSendVars, numProcs * sizeof(int));
290:     PetscMemzero(sumSendVars,   numProcs * sizeof(int));
291:     PetscMemzero(ghostRecvVars, numProcs * sizeof(int));
292:     PetscMemzero(sumRecvVars,   numProcs * sizeof(int));
293:     PetscMemzero(displs,        numProcs * sizeof(int));

295:     /* Get number of ghost variables to receive from each processor and size of blocks --
296:          we here assume that classes[] already has ghost node classes in it */
297:     for(node = 0; node < numGhostNodes; node++) {
298:       PartitionGhostToGlobalNodeIndex(part, node, &gNode, &proc);
299:       nclass = classes[numNodes+node];
300:       ghostRecvVars[proc]++;
301:       o->numOverlapVars += constMap->classSizes[nclass];
302:     }

304:     /* Get number of constrained ghost variables to receive from each processor and size of blocks */
305:     (*constCtx->ops->getsize)(constCtx, PETSC_NULL, PETSC_NULL, &o->numOverlapNewVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
306: 
307:     o->numOverlapVars += o->numOverlapNewVars - o->numLocNewVars;

309:     /* Get sizes of ghost variable blocks to send to each processor */
310:     MPI_Alltoall(ghostRecvVars, 1, MPI_INT, ghostSendVars, 1, MPI_INT, o->comm);

312:     /* Calculate offets into the ghost variable receive array */
313:     for(proc = 1; proc < numProcs; proc++) {
314:       sumRecvVars[proc] = sumRecvVars[proc-1] + ghostRecvVars[proc-1];
315:       displs[proc]      = sumRecvVars[proc];
316:     }

318:     /* Calculate offsets into the ghost variable send array */
319:     for(proc = 1; proc < numProcs; proc++)
320:       sumSendVars[proc] = sumSendVars[proc-1] + ghostSendVars[proc-1];

322:     /* Send requests for ghost variable offsets to each processor */
323:     numSendGhostVars = sumSendVars[numProcs-1] + ghostSendVars[numProcs-1];
324:     PetscMalloc(numSendGhostVars * sizeof(int), &sendGhostBuffer);
325:     for(node = 0; node < numGhostNodes; node++) {
326:       PartitionGhostToGlobalNodeIndex(part, node, &gNode, &proc);
327:       o->offsets[numNodes+(displs[proc]++)] = gNode;
328:     }
329:     MPI_Alltoallv(&o->offsets[numNodes],  ghostRecvVars, sumRecvVars, MPI_INT,
330:                          sendGhostBuffer,        ghostSendVars, sumSendVars, MPI_INT, o->comm);
331: 

333:     /* Send ghost variables offsets to each processor */
334:     for(node = 0; node < numSendGhostVars; node++) {
335:       PartitionGlobalToLocalNodeIndex(part, sendGhostBuffer[node], &locNode);
336:       sendGhostBuffer[node] = o->offsets[locNode];
337:     }
338:     MPI_Alltoallv(sendGhostBuffer,       ghostSendVars, sumSendVars, MPI_INT,
339:                          &o->offsets[numNodes], ghostRecvVars, sumRecvVars, MPI_INT, o->comm);
340: 

342:     /* Cleanup */
343:     PetscFree(ghostSendVars);
344:     PetscFree(sumSendVars);
345:     PetscFree(ghostRecvVars);
346:     PetscFree(sumRecvVars);
347:     PetscFree(displs);
348:     PetscFree(sendGhostBuffer);

350:     /* We maintain local offsets for ghost variables, meaning the offsets after the last
351:        interior variable, rather than the offset of the given ghost variable in the global
352:        matrix. */
353:     PetscMalloc(numGhostNodes * sizeof(int), &o->localOffsets);
354:     for(node = 0, var = o->numLocVars; node < numGhostNodes; node++) {
355:       o->localOffsets[node] = var;
356:       nclass = classes[numNodes+node];
357:       var   += constMap->classSizes[nclass];
358:     }
359:   }

361:   /* Setup local field offsets */
362:   VarOrderingGetClassMap(oldOrder, &map);
363:   numClassesOld = map->numClasses;
364:   for(f = 0; f < numFields; f++) {
365:     field = fields[f];
366:     PetscMalloc(numClasses * sizeof(int), &o->localStart[field]);
367:     for(nclass = 0; nclass < numClassesOld; nclass++) {
368:       o->localStart[field][nclass] = localStart[field][nclass];
369:     }
370:     for(i = numClassesOld; i < numClasses; i++) {
371:       /* Invert the class map */
372:       for(nclass = 0; nclass < numClassesOld; nclass++) {
373:         if (classMap[nclass] == i) break;
374:       }
375:       if (nclass >= numClassesOld) SETERRQ(PETSC_ERR_PLIB, "Invalid class map");

377:       /* Subtract out the constrained fields */
378:       o->localStart[field][i]    = localStart[field][nclass];
379:       if (localStart[field][nclass] > localStart[constField][nclass])
380:         o->localStart[field][i] -= comp;
381:     }
382:   }

384:   *order = o;
385:   return(0);
386: }

388: #undef  __FUNCT__
390: int GridCreateVarScatter_Triangular_2D(Grid grid, VarOrdering order, GVec vec, VarOrdering embedOrder, GVec embedVec,
391:                                        VecScatter *scatter)
392: {
393:   PetscConstraintObject constCtx = grid->constraintCtx;
394:   FieldClassMap         map, embedMap;
395:   int                   numEmbedFields, numNodes;
396:   int                  *embedFields, *classes, *embedClassSizes;
397:   int                   numLocVars, numEmbedLocVars, numNewLocVars;
398:   int                  *offsets, **localStart;
399:   IS                    is, embedIS;
400:   int                  *indices;
401:   PetscTruth            isConstrained;
402:   int                   node, nclass, f, field, comp, var, count;
403:   int                   ierr;

406:   /* Retrieve orderings */
407:   VarOrderingGetClassMap(order,      &map);
408:   VarOrderingGetClassMap(embedOrder, &embedMap);
409:   numNodes        = map->numNodes;
410:   classes         = map->classes;
411:   numLocVars      = order->numLocVars;
412:   offsets         = order->offsets;
413:   localStart      = order->localStart;
414:   numEmbedFields  = embedMap->numFields;
415:   embedFields     = embedMap->fields;
416:   embedClassSizes = embedMap->classSizes;
417:   numEmbedLocVars = embedOrder->numLocVars;

419:   PetscMalloc(numEmbedLocVars * sizeof(int), &indices);
420:   for(node = 0, count = 0; node < numNodes; node++) {
421:     nclass = classes[node];
422:     if (embedClassSizes[nclass] > 0) {
423:       for(f = 0; f < numEmbedFields; f++) {
424:         field = embedFields[f];
425:         if (localStart[field][nclass] >= 0) {
426:           GridGetFieldComponents(grid, field, &comp);
427:           for(var = 0; var < comp; var++)
428:             indices[count++] = offsets[node] + localStart[field][nclass] + var;
429:         }
430:       }
431:     }
432:   }
433:   /* Handle new variables */
434:   GridIsConstrained(grid, &isConstrained);
435:   if ((isConstrained == PETSC_TRUE) && (count < numEmbedLocVars)) {
436:     (*constCtx->ops->getsize)(constCtx, &numNewLocVars, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL);
437: 
438:     if (count != numEmbedLocVars - numNewLocVars) SETERRQ(PETSC_ERR_PLIB, "Inconsistent variable orderings");
439:     /* Fill in the embed variables */
440:     for(var = numNewLocVars; var > 0; var--)
441:       indices[count++] = numLocVars - var;
442:   }
443:   if (count != numEmbedLocVars) SETERRQ(PETSC_ERR_PLIB, "Inconsistent variable orderings");

445:   /* Create mappings */
446:   ISCreateGeneral(order->comm, numEmbedLocVars, indices, &is);
447:   ISCreateStride(order->comm, numEmbedLocVars, 0, 1, &embedIS);

449:   /* Create the scatter */
450:   VecScatterCreate(vec, is, embedVec, embedIS, scatter);

452:   /* Cleanup */
453:   ISDestroy(is);
454:   ISDestroy(embedIS);
455:   PetscFree(indices);

457:   return(0);
458: }

460: #undef  __FUNCT__
462: int GridCreateLocalVarOrdering_Triangular_2D(Grid grid, int numFields, int *fields, LocalVarOrdering *order)
463: {
464:   LocalVarOrdering o;
465:   int              elemOffset;
466:   int              f, field;
467:   int              ierr;

470:   /* Create the ordering */
471:   PetscHeaderCreate(o, _LocalVarOrdering, int, VAR_ORDER_COOKIE, 0, "LocalVarOrdering", grid->comm, LocalVarOrderingDestroy, 0);
472:   PetscLogObjectCreate(o);

474:   /* Allocate memory */
475:   o->numFields = numFields;
476:   PetscMalloc(numFields    * sizeof(int), &o->fields);
477:   PetscMalloc(grid->numFields * sizeof(int), &o->elemStart);
478:   PetscLogObjectMemory(o, (numFields + grid->numFields) * sizeof(int));
479:   PetscMemcpy(o->fields, fields, numFields * sizeof(int));

481:   /* Put in sentinel values */
482:   for(f = 0; f < grid->numFields; f++) {
483:     o->elemStart[f] = -1;
484:   }

486:   /* Setup local and global offsets offsets */
487:   for(f = 0, elemOffset = 0; f < numFields; f++) {
488:     field               = fields[f];
489:     o->elemStart[field] = elemOffset;
490:     elemOffset         += grid->fields[field].disc->size;
491:   }
492:   o->elemSize = elemOffset;
493:   *order = o;

495:   return(0);
496: }