xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 6ac5842e34eedc6428162d8d42bedaaf46eae34c)
1 
2 #include <petsc-private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
3 #include <petscdmda.h>
4 #include <petscviewer.h>
5 
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;
15 
16 PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
17 PetscInt       DMDAGetNeighborRelative(DM, PassiveReal, PassiveReal);
18 PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);
19 
20 PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
21 PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "CharacteristicView"
25 PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
26 {
27   PetscBool      iascii;
28   PetscErrorCode ierr;
29 
30   PetscFunctionBegin;
31   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
32   if (!viewer) {
33     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer);CHKERRQ(ierr);
34   }
35   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
36   PetscCheckSameComm(c, 1, viewer, 2);
37 
38   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
39   if (!iascii) {
40     if (c->ops->view) {
41       ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr);
42     }
43   }
44   PetscFunctionReturn(0);
45 }
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "CharacteristicDestroy"
49 PetscErrorCode CharacteristicDestroy(Characteristic *c)
50 {
51   PetscErrorCode ierr;
52 
53   PetscFunctionBegin;
54   if (!*c) PetscFunctionReturn(0);
55   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
56   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
57 
58   if ((*c)->ops->destroy) {
59     ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr);
60   }
61   ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr);
62   ierr = PetscFree((*c)->queue);CHKERRQ(ierr);
63   ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr);
64   ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr);
65   ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr);
66   ierr = PetscFree((*c)->needCount);CHKERRQ(ierr);
67   ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr);
68   ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr);
69   ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr);
70   ierr = PetscFree((*c)->request);CHKERRQ(ierr);
71   ierr = PetscFree((*c)->status);CHKERRQ(ierr);
72   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "CharacteristicCreate"
78 PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
79 {
80   Characteristic newC;
81   PetscErrorCode ierr;
82 
83   PetscFunctionBegin;
84   PetscValidPointer(c, 2);
85   *c = NULL;
86 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
87   ierr = CharacteristicInitializePackage(NULL);CHKERRQ(ierr);
88 #endif
89 
90   ierr = PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "SemiLagrange", comm, CharacteristicDestroy, CharacteristicView);CHKERRQ(ierr);
91   ierr = PetscLogObjectCreate(newC);CHKERRQ(ierr);
92   *c   = newC;
93 
94   newC->structured          = PETSC_TRUE;
95   newC->numIds              = 0;
96   newC->velocityDA          = NULL;
97   newC->velocity            = NULL;
98   newC->velocityOld         = NULL;
99   newC->numVelocityComp     = 0;
100   newC->velocityComp        = NULL;
101   newC->velocityInterp      = NULL;
102   newC->velocityInterpLocal = NULL;
103   newC->velocityCtx         = NULL;
104   newC->fieldDA             = NULL;
105   newC->field               = NULL;
106   newC->numFieldComp        = 0;
107   newC->fieldComp           = NULL;
108   newC->fieldInterp         = NULL;
109   newC->fieldInterpLocal    = NULL;
110   newC->fieldCtx            = NULL;
111   newC->itemType            = 0;
112   newC->queue               = NULL;
113   newC->queueSize           = 0;
114   newC->queueMax            = 0;
115   newC->queueLocal          = NULL;
116   newC->queueLocalSize      = 0;
117   newC->queueLocalMax       = 0;
118   newC->queueRemote         = NULL;
119   newC->queueRemoteSize     = 0;
120   newC->queueRemoteMax      = 0;
121   newC->numNeighbors        = 0;
122   newC->neighbors           = NULL;
123   newC->needCount           = NULL;
124   newC->localOffsets        = NULL;
125   newC->fillCount           = NULL;
126   newC->remoteOffsets       = NULL;
127   newC->request             = NULL;
128   newC->status              = NULL;
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "CharacteristicSetType"
134 /*@C
135    CharacteristicSetType - Builds Characteristic for a particular solver.
136 
137    Logically Collective on Characteristic
138 
139    Input Parameters:
140 +  c    - the method of characteristics context
141 -  type - a known method
142 
143    Options Database Key:
144 .  -characteristic_type <method> - Sets the method; use -help for a list
145     of available methods
146 
147    Notes:
148    See "include/petsccharacteristic.h" for available methods
149 
150   Normally, it is best to use the CharacteristicSetFromOptions() command and
151   then set the Characteristic type from the options database rather than by using
152   this routine.  Using the options database provides the user with
153   maximum flexibility in evaluating the many different Krylov methods.
154   The CharacteristicSetType() routine is provided for those situations where it
155   is necessary to set the iterative solver independently of the command
156   line or options database.  This might be the case, for example, when
157   the choice of iterative solver changes during the execution of the
158   program, and the user's application is taking responsibility for
159   choosing the appropriate method.  In other words, this routine is
160   not for beginners.
161 
162   Level: intermediate
163 
164 .keywords: Characteristic, set, method
165 
166 .seealso: CharacteristicType
167 
168 @*/
169 PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
170 {
171   PetscErrorCode ierr, (*r)(Characteristic);
172   PetscBool      match;
173 
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
176   PetscValidCharPointer(type, 2);
177 
178   ierr = PetscObjectTypeCompare((PetscObject) c, type, &match);CHKERRQ(ierr);
179   if (match) PetscFunctionReturn(0);
180 
181   if (c->data) {
182     /* destroy the old private Characteristic context */
183     ierr = (*c->ops->destroy)(c);CHKERRQ(ierr);
184     c->ops->destroy = NULL;
185     c->data         = 0;
186   }
187 
188   ierr =  PetscFunctionListFind(PetscObjectComm((PetscObject)c),CharacteristicList, type,PETSC_TRUE, (void (**)(void)) &r);CHKERRQ(ierr);
189   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
190   c->setupcalled = 0;
191   ierr = (*r)(c);CHKERRQ(ierr);
192   ierr = PetscObjectChangeTypeName((PetscObject) c, type);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "CharacteristicSetUp"
198 /*@
199    CharacteristicSetUp - Sets up the internal data structures for the
200    later use of an iterative solver.
201 
202    Collective on Characteristic
203 
204    Input Parameter:
205 .  ksp   - iterative context obtained from CharacteristicCreate()
206 
207    Level: developer
208 
209 .keywords: Characteristic, setup
210 
211 .seealso: CharacteristicCreate(), CharacteristicSolve(), CharacteristicDestroy()
212 @*/
213 PetscErrorCode CharacteristicSetUp(Characteristic c)
214 {
215   PetscErrorCode ierr;
216 
217   PetscFunctionBegin;
218   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
219 
220   if (!((PetscObject)c)->type_name) {
221     ierr = CharacteristicSetType(c, CHARACTERISTICDA);CHKERRQ(ierr);
222   }
223 
224   if (c->setupcalled == 2) PetscFunctionReturn(0);
225 
226   ierr = PetscLogEventBegin(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
227   if (!c->setupcalled) {
228     ierr = (*c->ops->setup)(c);CHKERRQ(ierr);
229   }
230   ierr = PetscLogEventEnd(CHARACTERISTIC_SetUp,c,0,0,0);CHKERRQ(ierr);
231   c->setupcalled = 2;
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "CharacteristicRegister"
237 /*@C
238   CharacteristicRegister - See CharacteristicRegisterDynamic()
239 
240   Level: advanced
241 @*/
242 PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic))
243 {
244   PetscErrorCode ierr;
245   char           fullname[PETSC_MAX_PATH_LEN];
246 
247   PetscFunctionBegin;
248   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
249   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "CharacteristicSetVelocityInterpolation"
255 PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
256 {
257   PetscFunctionBegin;
258   c->velocityDA      = da;
259   c->velocity        = v;
260   c->velocityOld     = vOld;
261   c->numVelocityComp = numComponents;
262   c->velocityComp    = components;
263   c->velocityInterp  = interp;
264   c->velocityCtx     = ctx;
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal"
270 PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
271 {
272   PetscFunctionBegin;
273   c->velocityDA          = da;
274   c->velocity            = v;
275   c->velocityOld         = vOld;
276   c->numVelocityComp     = numComponents;
277   c->velocityComp        = components;
278   c->velocityInterpLocal = interp;
279   c->velocityCtx         = ctx;
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "CharacteristicSetFieldInterpolation"
285 PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
286 {
287   PetscFunctionBegin;
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   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal"
302 PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
303 {
304   PetscFunctionBegin;
305 #if 0
306   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.");
307 #endif
308   c->fieldDA          = da;
309   c->field            = v;
310   c->numFieldComp     = numComponents;
311   c->fieldComp        = components;
312   c->fieldInterpLocal = interp;
313   c->fieldCtx         = ctx;
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "CharacteristicSolve"
319 PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
320 {
321   CharacteristicPointDA2D Qi;
322   DM                      da = c->velocityDA;
323   Vec                     velocityLocal, velocityLocalOld;
324   Vec                     fieldLocal;
325   DMDALocalInfo           info;
326   PetscScalar             **solArray;
327   void                    *velocityArray;
328   void                    *velocityArrayOld;
329   void                    *fieldArray;
330   PassiveScalar           *interpIndices;
331   PassiveScalar           *velocityValues, *velocityValuesOld;
332   PassiveScalar           *fieldValues;
333   PetscMPIInt             rank;
334   PetscInt                dim;
335   PetscMPIInt             neighbors[9];
336   PetscInt                dof;
337   PetscInt                gx, gy;
338   PetscInt                n, is, ie, js, je, comp;
339   PetscErrorCode          ierr;
340 
341   PetscFunctionBegin;
342   c->queueSize = 0;
343   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
344   ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr);
345   ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr);
346   ierr = CharacteristicSetUp(c);CHKERRQ(ierr);
347   /* global and local grid info */
348   ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr);
349   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
350   is   = info.xs;          ie   = info.xs+info.xm;
351   js   = info.ys;          je   = info.ys+info.ym;
352   /* Allocation */
353   ierr = PetscMalloc(dim*sizeof(PetscScalar),                &interpIndices);CHKERRQ(ierr);
354   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr);
355   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr);
356   ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar),    &fieldValues);CHKERRQ(ierr);
357   ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
358 
359   /* -----------------------------------------------------------------------
360      PART 1, AT t-dt/2
361      -----------------------------------------------------------------------*/
362   ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
363   /* GET POSITION AT HALF TIME IN THE PAST */
364   if (c->velocityInterpLocal) {
365     ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
366     ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
367     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
368     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
369     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
370     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
371     ierr = DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
372     ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
373   }
374   ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr);
375   for (Qi.j = js; Qi.j < je; Qi.j++) {
376     for (Qi.i = is; Qi.i < ie; Qi.i++) {
377       interpIndices[0] = Qi.i;
378       interpIndices[1] = Qi.j;
379       if (c->velocityInterpLocal) c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
380       else c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
381       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
382       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
383 
384       /* Determine whether the position at t - dt/2 is local */
385       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
386 
387       /* Check for Periodic boundaries and move all periodic points back onto the domain */
388       ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
389       ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr);
390     }
391   }
392   ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
393 
394   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
395   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
396   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
397 
398   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
399   /* Calculate velocity at t_n+1/2 (local values) */
400   ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr);
401   for (n = 0; n < c->queueSize; n++) {
402     Qi = c->queue[n];
403     if (c->neighbors[Qi.proc] == rank) {
404       interpIndices[0] = Qi.x;
405       interpIndices[1] = Qi.y;
406       if (c->velocityInterpLocal) {
407         c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
408         c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
409       } else {
410         c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
411         c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
412       }
413       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
414       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
415     }
416     c->queue[n] = Qi;
417   }
418   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
419 
420   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
421   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
422   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
423 
424 
425   /* Calculate velocity at t_n+1/2 (fill remote requests) */
426   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
427   ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr);
428   for (n = 0; n < c->queueRemoteSize; n++) {
429     Qi = c->queueRemote[n];
430     interpIndices[0] = Qi.x;
431     interpIndices[1] = Qi.y;
432     if (c->velocityInterpLocal) {
433       c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
434       c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
435     } else {
436       c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
437       c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
438     }
439     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
440     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
441     c->queueRemote[n] = Qi;
442   }
443   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
444   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
445   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
446   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
447   if (c->velocityInterpLocal) {
448     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
449     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
450     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
451     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
452   }
453   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
454 
455   /* -----------------------------------------------------------------------
456      PART 2, AT t-dt
457      -----------------------------------------------------------------------*/
458 
459   /* GET POSITION AT t_n (local values) */
460   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
461   ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr);
462   for (n = 0; n < c->queueSize; n++) {
463     Qi   = c->queue[n];
464     Qi.x = Qi.i - Qi.x*dt;
465     Qi.y = Qi.j - Qi.y*dt;
466 
467     /* Determine whether the position at t-dt is local */
468     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
469 
470     /* Check for Periodic boundaries and move all periodic points back onto the domain */
471     ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
472 
473     c->queue[n] = Qi;
474   }
475   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
476 
477   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
478   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
479   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
480 
481   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
482   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
483   if (c->fieldInterpLocal) {
484     ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
485     ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
486     ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
487     ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
488   }
489   ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr);
490   for (n = 0; n < c->queueSize; n++) {
491     if (c->neighbors[c->queue[n].proc] == rank) {
492       interpIndices[0] = c->queue[n].x;
493       interpIndices[1] = c->queue[n].y;
494       if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
495       else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
496       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
497     }
498   }
499   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
500 
501   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
502   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
503   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
504 
505   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
506   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
507   ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr);
508   for (n = 0; n < c->queueRemoteSize; n++) {
509     interpIndices[0] = c->queueRemote[n].x;
510     interpIndices[1] = c->queueRemote[n].y;
511 
512     /* for debugging purposes */
513     if (1) { /* hacked bounds test...let's do better */
514       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
515 
516       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);
517     }
518 
519     if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
520     else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
521     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
522   }
523   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
524 
525   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
526   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
527   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
528   if (c->fieldInterpLocal) {
529     ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
530     ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
531   }
532   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
533 
534   /* Return field of characteristics at t_n-1 */
535   ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
536   ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
537   ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
538   for (n = 0; n < c->queueSize; n++) {
539     Qi = c->queue[n];
540     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
541   }
542   ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
543   ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
544   ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
545 
546   /* Cleanup */
547   ierr = PetscFree(interpIndices);CHKERRQ(ierr);
548   ierr = PetscFree(velocityValues);CHKERRQ(ierr);
549   ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr);
550   ierr = PetscFree(fieldValues);CHKERRQ(ierr);
551   PetscFunctionReturn(0);
552 }
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "CharacteristicSetNeighbors"
556 PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
557 {
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   c->numNeighbors = numNeighbors;
562   ierr = PetscFree(c->neighbors);CHKERRQ(ierr);
563   ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr);
564   ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr);
565   PetscFunctionReturn(0);
566 }
567 
568 #undef __FUNCT__
569 #define __FUNCT__ "CharacteristicAddPoint"
570 PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
571 {
572   PetscFunctionBegin;
573   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
574   c->queue[c->queueSize++] = *point;
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "CharacteristicSendCoordinatesBegin"
580 int CharacteristicSendCoordinatesBegin(Characteristic c)
581 {
582   PetscMPIInt    rank, tag = 121;
583   PetscInt       i, n;
584   PetscErrorCode ierr;
585 
586   PetscFunctionBegin;
587   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
588   ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr);
589   ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr);
590   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
591   c->fillCount[0] = 0;
592   for (n = 1; n < c->numNeighbors; n++) {
593     ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
594   }
595   for (n = 1; n < c->numNeighbors; n++) {
596     ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
597   }
598   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
599   /* Initialize the remote queue */
600   c->queueLocalMax  = c->localOffsets[0]  = 0;
601   c->queueRemoteMax = c->remoteOffsets[0] = 0;
602   for (n = 1; n < c->numNeighbors; n++) {
603     c->remoteOffsets[n] = c->queueRemoteMax;
604     c->queueRemoteMax  += c->fillCount[n];
605     c->localOffsets[n]  = c->queueLocalMax;
606     c->queueLocalMax   += c->needCount[n];
607   }
608   /* HACK BEGIN */
609   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
610   c->needCount[0] = 0;
611   /* HACK END */
612   if (c->queueRemoteMax) {
613     ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr);
614   } else c->queueRemote = NULL;
615   c->queueRemoteSize = c->queueRemoteMax;
616 
617   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
618   for (n = 1; n < c->numNeighbors; n++) {
619     ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr);
620     ierr = MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
621   }
622   for (n = 1; n < c->numNeighbors; n++) {
623     ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr);
624     ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
625   }
626   PetscFunctionReturn(0);
627 }
628 
629 #undef __FUNCT__
630 #define __FUNCT__ "CharacteristicSendCoordinatesEnd"
631 PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
632 {
633 #if 0
634   PetscMPIInt rank;
635   PetscInt    n;
636 #endif
637   PetscErrorCode ierr;
638 
639   PetscFunctionBegin;
640   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
641 #if 0
642   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
643   for (n = 0; n < c->queueRemoteSize; n++) {
644     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);
645   }
646 #endif
647   PetscFunctionReturn(0);
648 }
649 
650 #undef __FUNCT__
651 #define __FUNCT__ "CharacteristicGetValuesBegin"
652 PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
653 {
654   PetscMPIInt    tag = 121;
655   PetscInt       n;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
660   for (n = 1; n < c->numNeighbors; n++) {
661     ierr = MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
662   }
663   for (n = 1; n < c->numNeighbors; n++) {
664     ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
665   }
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "CharacteristicGetValuesEnd"
671 PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
672 {
673   PetscErrorCode ierr;
674 
675   PetscFunctionBegin;
676   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
677   /* Free queue of requests from other procs */
678   ierr = PetscFree(c->queueRemote);CHKERRQ(ierr);
679   PetscFunctionReturn(0);
680 }
681 
682 /*---------------------------------------------------------------------*/
683 #undef __FUNCT__
684 #define __FUNCT__ "CharacteristicHeapSort"
685 /*
686   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
687 */
688 PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
689 /*---------------------------------------------------------------------*/
690 {
691   PetscErrorCode          ierr;
692   CharacteristicPointDA2D temp;
693   PetscInt                n;
694 
695   PetscFunctionBegin;
696   if (0) { /* Check the order of the queue before sorting */
697     PetscInfo(NULL, "Before Heap sort\n");
698     for (n=0; n<size; n++) {
699       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
700     }
701   }
702 
703   /* SORTING PHASE */
704   for (n = (size / 2)-1; n >= 0; n--) {
705     ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/
706   }
707   for (n = size-1; n >= 1; n--) {
708     temp     = queue[0];
709     queue[0] = queue[n];
710     queue[n] = temp;
711     ierr     = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr);
712   }
713   if (0) { /* Check the order of the queue after sorting */
714     ierr = PetscInfo(NULL, "Avter  Heap sort\n");CHKERRQ(ierr);
715     for (n=0; n<size; n++) {
716       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
717     }
718   }
719   PetscFunctionReturn(0);
720 }
721 
722 /*---------------------------------------------------------------------*/
723 #undef __FUNCT__
724 #define __FUNCT__ "CharacteristicSiftDown"
725 /*
726   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
727 */
728 PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
729 /*---------------------------------------------------------------------*/
730 {
731   PetscBool               done = PETSC_FALSE;
732   PetscInt                maxChild;
733   CharacteristicPointDA2D temp;
734 
735   PetscFunctionBegin;
736   while ((root*2 <= bottom) && (!done)) {
737     if (root*2 == bottom) maxChild = root * 2;
738     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
739     else maxChild = root * 2 + 1;
740 
741     if (queue[root].proc < queue[maxChild].proc) {
742       temp            = queue[root];
743       queue[root]     = queue[maxChild];
744       queue[maxChild] = temp;
745       root            = maxChild;
746     } else done = PETSC_TRUE;
747   }
748   PetscFunctionReturn(0);
749 }
750 
751 #undef __FUNCT__
752 #define __FUNCT__ "DMDAGetNeighborsRank"
753 /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
754 PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
755 {
756   DMDABoundaryType bx, by;
757   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
758   MPI_Comm         comm;
759   PetscMPIInt      rank;
760   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
761   PetscErrorCode   ierr;
762 
763   PetscFunctionBegin;
764   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
765   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
766   ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr);
767 
768   if (bx == DMDA_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
769   if (by == DMDA_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
770 
771   neighbors[0] = rank;
772   rank = 0;
773   ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr);
774   for (pj=0; pj<PJ; pj++) {
775     ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr);
776     for (pi=0; pi<PI; pi++) {
777       procs[pj][pi] = rank;
778       rank++;
779     }
780   }
781 
782   pi  = neighbors[0] % PI;
783   pj  = neighbors[0] / PI;
784   pim = pi-1;  if (pim<0) pim=PI-1;
785   pip = (pi+1)%PI;
786   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
787   pjp = (pj+1)%PJ;
788 
789   neighbors[1] = procs[pj] [pim];
790   neighbors[2] = procs[pjp][pim];
791   neighbors[3] = procs[pjp][pi];
792   neighbors[4] = procs[pjp][pip];
793   neighbors[5] = procs[pj] [pip];
794   neighbors[6] = procs[pjm][pip];
795   neighbors[7] = procs[pjm][pi];
796   neighbors[8] = procs[pjm][pim];
797 
798   if (!IPeriodic) {
799     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
800     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
801   }
802 
803   if (!JPeriodic) {
804     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
805     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
806   }
807 
808   for (pj = 0; pj < PJ; pj++) {
809     ierr = PetscFree(procs[pj]);CHKERRQ(ierr);
810   }
811   ierr = PetscFree(procs);CHKERRQ(ierr);
812   PetscFunctionReturn(0);
813 }
814 
815 #undef __FUNCT__
816 #define __FUNCT__ "DMDAGetNeighborRelative"
817 /*
818   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
819     2 | 3 | 4
820     __|___|__
821     1 | 0 | 5
822     __|___|__
823     8 | 7 | 6
824       |   |
825 */
826 PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr)
827 {
828   DMDALocalInfo  info;
829   PassiveReal    is,ie,js,je;
830   PetscErrorCode ierr;
831 
832   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
833   is   = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5;
834   js   = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5;
835 
836   if (ir >= is && ir <= ie) { /* center column */
837     if (jr >= js && jr <= je) return 0;
838     else if (jr < js)         return 7;
839     else                      return 3;
840   } else if (ir < is) {     /* left column */
841     if (jr >= js && jr <= je) return 1;
842     else if (jr < js)         return 8;
843     else                      return 2;
844   } else {                  /* right column */
845     if (jr >= js && jr <= je) return 5;
846     else if (jr < js)         return 6;
847     else                      return 4;
848   }
849 }
850