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