Actual source code: characteristic.c


  2: #include <petsc/private/characteristicimpl.h>
  3: #include <petscdmda.h>
  4: #include <petscviewer.h>

  6: PetscClassId  CHARACTERISTIC_CLASSID;
  7: PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
  8: PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
  9: PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
 10: /*
 11:    Contains the list of registered characteristic routines
 12: */
 13: PetscFunctionList CharacteristicList              = NULL;
 14: PetscBool         CharacteristicRegisterAllCalled = PETSC_FALSE;

 16: PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
 17: PetscInt       DMDAGetNeighborRelative(DM, PetscReal, PetscReal);
 18: PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);

 20: PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
 21: PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);

 23: PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
 24: {
 25:   PetscBool      iascii;

 30:   if (!viewer) {
 31:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);
 32:   }

 36:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 37:   if (!iascii) {
 38:     if (c->ops->view) {
 39:       (*c->ops->view)(c, viewer);
 40:     }
 41:   }
 42:   return(0);
 43: }

 45: PetscErrorCode CharacteristicDestroy(Characteristic *c)
 46: {

 50:   if (!*c) return(0);
 52:   if (--((PetscObject)(*c))->refct > 0) return(0);

 54:   if ((*c)->ops->destroy) {
 55:     (*(*c)->ops->destroy)((*c));
 56:   }
 57:   MPI_Type_free(&(*c)->itemType);
 58:   PetscFree((*c)->queue);
 59:   PetscFree((*c)->queueLocal);
 60:   PetscFree((*c)->queueRemote);
 61:   PetscFree((*c)->neighbors);
 62:   PetscFree((*c)->needCount);
 63:   PetscFree((*c)->localOffsets);
 64:   PetscFree((*c)->fillCount);
 65:   PetscFree((*c)->remoteOffsets);
 66:   PetscFree((*c)->request);
 67:   PetscFree((*c)->status);
 68:   PetscHeaderDestroy(c);
 69:   return(0);
 70: }

 72: PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
 73: {
 74:   Characteristic newC;

 79:   *c = NULL;
 80:   CharacteristicInitializePackage();

 82:   PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView);
 83:   *c   = newC;

 85:   newC->structured          = PETSC_TRUE;
 86:   newC->numIds              = 0;
 87:   newC->velocityDA          = NULL;
 88:   newC->velocity            = NULL;
 89:   newC->velocityOld         = NULL;
 90:   newC->numVelocityComp     = 0;
 91:   newC->velocityComp        = NULL;
 92:   newC->velocityInterp      = NULL;
 93:   newC->velocityInterpLocal = NULL;
 94:   newC->velocityCtx         = NULL;
 95:   newC->fieldDA             = NULL;
 96:   newC->field               = NULL;
 97:   newC->numFieldComp        = 0;
 98:   newC->fieldComp           = NULL;
 99:   newC->fieldInterp         = NULL;
100:   newC->fieldInterpLocal    = NULL;
101:   newC->fieldCtx            = NULL;
102:   newC->itemType            = 0;
103:   newC->queue               = NULL;
104:   newC->queueSize           = 0;
105:   newC->queueMax            = 0;
106:   newC->queueLocal          = NULL;
107:   newC->queueLocalSize      = 0;
108:   newC->queueLocalMax       = 0;
109:   newC->queueRemote         = NULL;
110:   newC->queueRemoteSize     = 0;
111:   newC->queueRemoteMax      = 0;
112:   newC->numNeighbors        = 0;
113:   newC->neighbors           = NULL;
114:   newC->needCount           = NULL;
115:   newC->localOffsets        = NULL;
116:   newC->fillCount           = NULL;
117:   newC->remoteOffsets       = NULL;
118:   newC->request             = NULL;
119:   newC->status              = NULL;
120:   return(0);
121: }

123: /*@C
124:    CharacteristicSetType - Builds Characteristic for a particular solver.

126:    Logically Collective on Characteristic

128:    Input Parameters:
129: +  c    - the method of characteristics context
130: -  type - a known method

132:    Options Database Key:
133: .  -characteristic_type <method> - Sets the method; use -help for a list
134:     of available methods

136:    Notes:
137:    See "include/petsccharacteristic.h" for available methods

139:   Normally, it is best to use the CharacteristicSetFromOptions() command and
140:   then set the Characteristic type from the options database rather than by using
141:   this routine.  Using the options database provides the user with
142:   maximum flexibility in evaluating the many different Krylov methods.
143:   The CharacteristicSetType() routine is provided for those situations where it
144:   is necessary to set the iterative solver independently of the command
145:   line or options database.  This might be the case, for example, when
146:   the choice of iterative solver changes during the execution of the
147:   program, and the user's application is taking responsibility for
148:   choosing the appropriate method.  In other words, this routine is
149:   not for beginners.

151:   Level: intermediate

153: .seealso: CharacteristicType

155: @*/
156: PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
157: {
158:   PetscErrorCode ierr, (*r)(Characteristic);
159:   PetscBool      match;


165:   PetscObjectTypeCompare((PetscObject) c, type, &match);
166:   if (match) return(0);

168:   if (c->data) {
169:     /* destroy the old private Characteristic context */
170:     (*c->ops->destroy)(c);
171:     c->ops->destroy = NULL;
172:     c->data         = NULL;
173:   }

175:    PetscFunctionListFind(CharacteristicList,type,&r);
176:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
177:   c->setupcalled = 0;
178:   (*r)(c);
179:   PetscObjectChangeTypeName((PetscObject) c, type);
180:   return(0);
181: }

183: /*@
184:    CharacteristicSetUp - Sets up the internal data structures for the
185:    later use of an iterative solver.

187:    Collective on Characteristic

189:    Input Parameter:
190: .  ksp   - iterative context obtained from CharacteristicCreate()

192:    Level: developer

194: .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
195: @*/
196: PetscErrorCode CharacteristicSetUp(Characteristic c)
197: {


203:   if (!((PetscObject)c)->type_name) {
204:     CharacteristicSetType(c, CHARACTERISTICDA);
205:   }

207:   if (c->setupcalled == 2) return(0);

209:   PetscLogEventBegin(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL);
210:   if (!c->setupcalled) {
211:     (*c->ops->setup)(c);
212:   }
213:   PetscLogEventEnd(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL);
214:   c->setupcalled = 2;
215:   return(0);
216: }

218: /*@C
219:   CharacteristicRegister -  Adds a solver to the method of characteristics package.

221:    Not Collective

223:    Input Parameters:
224: +  name_solver - name of a new user-defined solver
225: -  routine_create - routine to create method context

227:   Sample usage:
228: .vb
229:     CharacteristicRegister("my_char", MyCharCreate);
230: .ve

232:   Then, your Characteristic type can be chosen with the procedural interface via
233: .vb
234:     CharacteristicCreate(MPI_Comm, Characteristic* &char);
235:     CharacteristicSetType(char,"my_char");
236: .ve
237:    or at runtime via the option
238: .vb
239:     -characteristic_type my_char
240: .ve

242:    Notes:
243:    CharacteristicRegister() may be called multiple times to add several user-defined solvers.

245: .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy()

247:   Level: advanced
248: @*/
249: PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic))
250: {

254:   CharacteristicInitializePackage();
255:   PetscFunctionListAdd(&CharacteristicList,sname,function);
256:   return(0);
257: }

259: PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
260: {
262:   c->velocityDA      = da;
263:   c->velocity        = v;
264:   c->velocityOld     = vOld;
265:   c->numVelocityComp = numComponents;
266:   c->velocityComp    = components;
267:   c->velocityInterp  = interp;
268:   c->velocityCtx     = ctx;
269:   return(0);
270: }

272: PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
273: {
275:   c->velocityDA          = da;
276:   c->velocity            = v;
277:   c->velocityOld         = vOld;
278:   c->numVelocityComp     = numComponents;
279:   c->velocityComp        = components;
280:   c->velocityInterpLocal = interp;
281:   c->velocityCtx         = ctx;
282:   return(0);
283: }

285: PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
286: {
288: #if 0
289:   if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
290: #endif
291:   c->fieldDA      = da;
292:   c->field        = v;
293:   c->numFieldComp = numComponents;
294:   c->fieldComp    = components;
295:   c->fieldInterp  = interp;
296:   c->fieldCtx     = ctx;
297:   return(0);
298: }

300: PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
301: {
303: #if 0
304:   if (numComponents > 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
305: #endif
306:   c->fieldDA          = da;
307:   c->field            = v;
308:   c->numFieldComp     = numComponents;
309:   c->fieldComp        = components;
310:   c->fieldInterpLocal = interp;
311:   c->fieldCtx         = ctx;
312:   return(0);
313: }

315: PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
316: {
317:   CharacteristicPointDA2D Qi;
318:   DM                      da = c->velocityDA;
319:   Vec                     velocityLocal, velocityLocalOld;
320:   Vec                     fieldLocal;
321:   DMDALocalInfo           info;
322:   PetscScalar             **solArray;
323:   void                    *velocityArray;
324:   void                    *velocityArrayOld;
325:   void                    *fieldArray;
326:   PetscScalar             *interpIndices;
327:   PetscScalar             *velocityValues, *velocityValuesOld;
328:   PetscScalar             *fieldValues;
329:   PetscMPIInt             rank;
330:   PetscInt                dim;
331:   PetscMPIInt             neighbors[9];
332:   PetscInt                dof;
333:   PetscInt                gx, gy;
334:   PetscInt                n, is, ie, js, je, comp;
335:   PetscErrorCode          ierr;

338:   c->queueSize = 0;
339:   MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);
340:   DMDAGetNeighborsRank(da, neighbors);
341:   CharacteristicSetNeighbors(c, 9, neighbors);
342:   CharacteristicSetUp(c);
343:   /* global and local grid info */
344:   DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
345:   DMDAGetLocalInfo(da, &info);
346:   is   = info.xs;          ie   = info.xs+info.xm;
347:   js   = info.ys;          je   = info.ys+info.ym;
348:   /* Allocation */
349:   PetscMalloc1(dim,                &interpIndices);
350:   PetscMalloc1(c->numVelocityComp, &velocityValues);
351:   PetscMalloc1(c->numVelocityComp, &velocityValuesOld);
352:   PetscMalloc1(c->numFieldComp,    &fieldValues);
353:   PetscLogEventBegin(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL);

355:   /* -----------------------------------------------------------------------
356:      PART 1, AT t-dt/2
357:      -----------------------------------------------------------------------*/
358:   PetscLogEventBegin(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL);
359:   /* GET POSITION AT HALF TIME IN THE PAST */
360:   if (c->velocityInterpLocal) {
361:     DMGetLocalVector(c->velocityDA, &velocityLocal);
362:     DMGetLocalVector(c->velocityDA, &velocityLocalOld);
363:     DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);
364:     DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);
365:     DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);
366:     DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);
367:     DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);
368:     DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);
369:   }
370:   PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");
371:   for (Qi.j = js; Qi.j < je; Qi.j++) {
372:     for (Qi.i = is; Qi.i < ie; Qi.i++) {
373:       interpIndices[0] = Qi.i;
374:       interpIndices[1] = Qi.j;
375:       if (c->velocityInterpLocal) {c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);}
376:       else {c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);}
377:       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
378:       Qi.y = Qi.j - velocityValues[1]*dt/2.0;

380:       /* Determine whether the position at t - dt/2 is local */
381:       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);

383:       /* Check for Periodic boundaries and move all periodic points back onto the domain */
384:       DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));
385:       CharacteristicAddPoint(c, &Qi);
386:     }
387:   }
388:   PetscLogEventEnd(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL);

390:   PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);
391:   CharacteristicSendCoordinatesBegin(c);
392:   PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);

394:   PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL);
395:   /* Calculate velocity at t_n+1/2 (local values) */
396:   PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");
397:   for (n = 0; n < c->queueSize; n++) {
398:     Qi = c->queue[n];
399:     if (c->neighbors[Qi.proc] == rank) {
400:       interpIndices[0] = Qi.x;
401:       interpIndices[1] = Qi.y;
402:       if (c->velocityInterpLocal) {
403:         c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
404:         c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
405:       } else {
406:         c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
407:         c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
408:       }
409:       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
410:       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
411:     }
412:     c->queue[n] = Qi;
413:   }
414:   PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL);

416:   PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);
417:   CharacteristicSendCoordinatesEnd(c);
418:   PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);

420:   /* Calculate velocity at t_n+1/2 (fill remote requests) */
421:   PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL);
422:   PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);
423:   for (n = 0; n < c->queueRemoteSize; n++) {
424:     Qi = c->queueRemote[n];
425:     interpIndices[0] = Qi.x;
426:     interpIndices[1] = Qi.y;
427:     if (c->velocityInterpLocal) {
428:       c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
429:       c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
430:     } else {
431:       c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
432:       c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
433:     }
434:     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
435:     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
436:     c->queueRemote[n] = Qi;
437:   }
438:   PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL);
439:   PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);
440:   CharacteristicGetValuesBegin(c);
441:   CharacteristicGetValuesEnd(c);
442:   if (c->velocityInterpLocal) {
443:     DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);
444:     DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);
445:     DMRestoreLocalVector(c->velocityDA, &velocityLocal);
446:     DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);
447:   }
448:   PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL);

450:   /* -----------------------------------------------------------------------
451:      PART 2, AT t-dt
452:      -----------------------------------------------------------------------*/

454:   /* GET POSITION AT t_n (local values) */
455:   PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL);
456:   PetscInfo(NULL, "Calculating position at t_{n}\n");
457:   for (n = 0; n < c->queueSize; n++) {
458:     Qi   = c->queue[n];
459:     Qi.x = Qi.i - Qi.x*dt;
460:     Qi.y = Qi.j - Qi.y*dt;

462:     /* Determine whether the position at t-dt is local */
463:     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);

465:     /* Check for Periodic boundaries and move all periodic points back onto the domain */
466:     DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));

468:     c->queue[n] = Qi;
469:   }
470:   PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL);

472:   PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);
473:   CharacteristicSendCoordinatesBegin(c);
474:   PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);

476:   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
477:   PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL);
478:   if (c->fieldInterpLocal) {
479:     DMGetLocalVector(c->fieldDA, &fieldLocal);
480:     DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);
481:     DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);
482:     DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);
483:   }
484:   PetscInfo(NULL, "Calculating local field at t_{n}\n");
485:   for (n = 0; n < c->queueSize; n++) {
486:     if (c->neighbors[c->queue[n].proc] == rank) {
487:       interpIndices[0] = c->queue[n].x;
488:       interpIndices[1] = c->queue[n].y;
489:       if (c->fieldInterpLocal) {c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);}
490:       else {c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);}
491:       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
492:     }
493:   }
494:   PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL);

496:   PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);
497:   CharacteristicSendCoordinatesEnd(c);
498:   PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);

500:   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
501:   PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL);
502:   PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);
503:   for (n = 0; n < c->queueRemoteSize; n++) {
504:     interpIndices[0] = c->queueRemote[n].x;
505:     interpIndices[1] = c->queueRemote[n].y;

507:     /* for debugging purposes */
508:     if (1) { /* hacked bounds test...let's do better */
509:       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];

511:       if ((im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar)  js - 1.) || (jm > (PetscScalar) je)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm);
512:     }

514:     if (c->fieldInterpLocal) {c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);}
515:     else {c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);}
516:     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
517:   }
518:   PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL);

520:   PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);
521:   CharacteristicGetValuesBegin(c);
522:   CharacteristicGetValuesEnd(c);
523:   if (c->fieldInterpLocal) {
524:     DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);
525:     DMRestoreLocalVector(c->fieldDA, &fieldLocal);
526:   }
527:   PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL);

529:   /* Return field of characteristics at t_n-1 */
530:   PetscLogEventBegin(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL);
531:   DMDAGetInfo(c->fieldDA,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
532:   DMDAVecGetArray(c->fieldDA, solution, &solArray);
533:   for (n = 0; n < c->queueSize; n++) {
534:     Qi = c->queue[n];
535:     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
536:   }
537:   DMDAVecRestoreArray(c->fieldDA, solution, &solArray);
538:   PetscLogEventEnd(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL);
539:   PetscLogEventEnd(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL);

541:   /* Cleanup */
542:   PetscFree(interpIndices);
543:   PetscFree(velocityValues);
544:   PetscFree(velocityValuesOld);
545:   PetscFree(fieldValues);
546:   return(0);
547: }

549: PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
550: {

554:   c->numNeighbors = numNeighbors;
555:   PetscFree(c->neighbors);
556:   PetscMalloc1(numNeighbors, &c->neighbors);
557:   PetscArraycpy(c->neighbors, neighbors, numNeighbors);
558:   return(0);
559: }

561: PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
562: {
564:   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %d", c->queueMax);
565:   c->queue[c->queueSize++] = *point;
566:   return(0);
567: }

569: int CharacteristicSendCoordinatesBegin(Characteristic c)
570: {
571:   PetscMPIInt    rank, tag = 121;
572:   PetscInt       i, n;

576:   MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);
577:   CharacteristicHeapSort(c, c->queue, c->queueSize);
578:   PetscArrayzero(c->needCount, c->numNeighbors);
579:   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
580:   c->fillCount[0] = 0;
581:   for (n = 1; n < c->numNeighbors; n++) {
582:     MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));
583:   }
584:   for (n = 1; n < c->numNeighbors; n++) {
585:     MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));
586:   }
587:   MPI_Waitall(c->numNeighbors-1, c->request, c->status);
588:   /* Initialize the remote queue */
589:   c->queueLocalMax  = c->localOffsets[0]  = 0;
590:   c->queueRemoteMax = c->remoteOffsets[0] = 0;
591:   for (n = 1; n < c->numNeighbors; n++) {
592:     c->remoteOffsets[n] = c->queueRemoteMax;
593:     c->queueRemoteMax  += c->fillCount[n];
594:     c->localOffsets[n]  = c->queueLocalMax;
595:     c->queueLocalMax   += c->needCount[n];
596:   }
597:   /* HACK BEGIN */
598:   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
599:   c->needCount[0] = 0;
600:   /* HACK END */
601:   if (c->queueRemoteMax) {
602:     PetscMalloc1(c->queueRemoteMax, &c->queueRemote);
603:   } else c->queueRemote = NULL;
604:   c->queueRemoteSize = c->queueRemoteMax;

606:   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
607:   for (n = 1; n < c->numNeighbors; n++) {
608:     PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);
609:     MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));
610:   }
611:   for (n = 1; n < c->numNeighbors; n++) {
612:     PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);
613:     MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));
614:   }
615:   return(0);
616: }

618: PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
619: {
620: #if 0
621:   PetscMPIInt rank;
622:   PetscInt    n;
623: #endif

627:   MPI_Waitall(c->numNeighbors-1, c->request, c->status);
628: #if 0
629:   MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);
630:   for (n = 0; n < c->queueRemoteSize; n++) {
631:     if (c->neighbors[c->queueRemote[n].proc] == rank) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is messed up, n = %d proc = %d", n, c->queueRemote[n].proc);
632:   }
633: #endif
634:   return(0);
635: }

637: PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
638: {
639:   PetscMPIInt    tag = 121;
640:   PetscInt       n;

644:   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
645:   for (n = 1; n < c->numNeighbors; n++) {
646:     MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));
647:   }
648:   for (n = 1; n < c->numNeighbors; n++) {
649:     MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));
650:   }
651:   return(0);
652: }

654: PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
655: {

659:   MPI_Waitall(c->numNeighbors-1, c->request, c->status);
660:   /* Free queue of requests from other procs */
661:   PetscFree(c->queueRemote);
662:   return(0);
663: }

665: /*---------------------------------------------------------------------*/
666: /*
667:   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
668: */
669: PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
670: /*---------------------------------------------------------------------*/
671: {
672:   PetscErrorCode          ierr;
673:   CharacteristicPointDA2D temp;
674:   PetscInt                n;

677:   if (0) { /* Check the order of the queue before sorting */
678:     PetscInfo(NULL, "Before Heap sort\n");
679:     for (n=0; n<size; n++) {
680:       PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);
681:     }
682:   }

684:   /* SORTING PHASE */
685:   for (n = (size / 2)-1; n >= 0; n--) {
686:     CharacteristicSiftDown(c, queue, n, size-1); /* Rich had size-1 here, Matt had size*/
687:   }
688:   for (n = size-1; n >= 1; n--) {
689:     temp     = queue[0];
690:     queue[0] = queue[n];
691:     queue[n] = temp;
692:     CharacteristicSiftDown(c, queue, 0, n-1);
693:   }
694:   if (0) { /* Check the order of the queue after sorting */
695:     PetscInfo(NULL, "Avter  Heap sort\n");
696:     for (n=0; n<size; n++) {
697:       PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);
698:     }
699:   }
700:   return(0);
701: }

703: /*---------------------------------------------------------------------*/
704: /*
705:   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
706: */
707: PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
708: /*---------------------------------------------------------------------*/
709: {
710:   PetscBool               done = PETSC_FALSE;
711:   PetscInt                maxChild;
712:   CharacteristicPointDA2D temp;

715:   while ((root*2 <= bottom) && (!done)) {
716:     if (root*2 == bottom) maxChild = root * 2;
717:     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
718:     else maxChild = root * 2 + 1;

720:     if (queue[root].proc < queue[maxChild].proc) {
721:       temp            = queue[root];
722:       queue[root]     = queue[maxChild];
723:       queue[maxChild] = temp;
724:       root            = maxChild;
725:     } else done = PETSC_TRUE;
726:   }
727:   return(0);
728: }

730: /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
731: PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
732: {
733:   DMBoundaryType   bx, by;
734:   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
735:   MPI_Comm         comm;
736:   PetscMPIInt      rank;
737:   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
738:   PetscErrorCode   ierr;

741:   PetscObjectGetComm((PetscObject) da, &comm);
742:   MPI_Comm_rank(comm, &rank);
743:   DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI,&PJ, NULL, NULL, NULL, &bx, &by,NULL, NULL);

745:   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
746:   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;

748:   neighbors[0] = rank;
749:   rank = 0;
750:   PetscMalloc1(PJ,&procs);
751:   for (pj=0; pj<PJ; pj++) {
752:     PetscMalloc1(PI,&(procs[pj]));
753:     for (pi=0; pi<PI; pi++) {
754:       procs[pj][pi] = rank;
755:       rank++;
756:     }
757:   }

759:   pi  = neighbors[0] % PI;
760:   pj  = neighbors[0] / PI;
761:   pim = pi-1;  if (pim<0) pim=PI-1;
762:   pip = (pi+1)%PI;
763:   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
764:   pjp = (pj+1)%PJ;

766:   neighbors[1] = procs[pj] [pim];
767:   neighbors[2] = procs[pjp][pim];
768:   neighbors[3] = procs[pjp][pi];
769:   neighbors[4] = procs[pjp][pip];
770:   neighbors[5] = procs[pj] [pip];
771:   neighbors[6] = procs[pjm][pip];
772:   neighbors[7] = procs[pjm][pi];
773:   neighbors[8] = procs[pjm][pim];

775:   if (!IPeriodic) {
776:     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
777:     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
778:   }

780:   if (!JPeriodic) {
781:     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
782:     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
783:   }

785:   for (pj = 0; pj < PJ; pj++) {
786:     PetscFree(procs[pj]);
787:   }
788:   PetscFree(procs);
789:   return(0);
790: }

792: /*
793:   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
794:     2 | 3 | 4
795:     __|___|__
796:     1 | 0 | 5
797:     __|___|__
798:     8 | 7 | 6
799:       |   |
800: */
801: PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
802: {
803:   DMDALocalInfo  info;
804:   PetscReal      is,ie,js,je;

807:   DMDAGetLocalInfo(da, &info);
808:   is   = (PetscReal) info.xs - 0.5; ie = (PetscReal) info.xs + info.xm - 0.5;
809:   js   = (PetscReal) info.ys - 0.5; je = (PetscReal) info.ys + info.ym - 0.5;

811:   if (ir >= is && ir <= ie) { /* center column */
812:     if (jr >= js && jr <= je) return 0;
813:     else if (jr < js)         return 7;
814:     else                      return 3;
815:   } else if (ir < is) {     /* left column */
816:     if (jr >= js && jr <= je) return 1;
817:     else if (jr < js)         return 8;
818:     else                      return 2;
819:   } else {                  /* right column */
820:     if (jr >= js && jr <= je) return 5;
821:     else if (jr < js)         return 6;
822:     else                      return 4;
823:   }
824: }