xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 009bbdc485cd9ad46be9940d3549e2dde9cdc322)
1 
2 #include <petsc-private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
3 
4 PetscClassId   CHARACTERISTIC_CLASSID;
5 PetscLogEvent  CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
6 PetscLogEvent  CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
7 PetscLogEvent  CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
8 /*
9    Contains the list of registered characteristic routines
10 */
11 PetscFunctionList  CharacteristicList = PETSC_NULL;
12 PetscBool          CharacteristicRegisterAllCalled = PETSC_FALSE;
13 
14 PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
15 PetscInt       DMDAGetNeighborRelative(DM, PassiveReal, PassiveReal);
16 PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);
17 
18 PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
19 PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "CharacteristicView"
23 PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24 {
25   PetscBool      iascii;
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
30   if (!viewer) {
31     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
32   }
33   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
34   PetscCheckSameComm(c, 1, viewer, 2);
35 
36   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
37   if (!iascii) {
38     if (c->ops->view) {
39       ierr = (*c->ops->view)(c, viewer);CHKERRQ(ierr);
40     }
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "CharacteristicDestroy"
47 PetscErrorCode CharacteristicDestroy(Characteristic *c)
48 {
49   PetscErrorCode ierr;
50 
51   PetscFunctionBegin;
52   if (!*c) PetscFunctionReturn(0);
53   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
54   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
55 
56   if ((*c)->ops->destroy) {
57     ierr = (*(*c)->ops->destroy)((*c));CHKERRQ(ierr);
58   }
59   ierr = MPI_Type_free(&(*c)->itemType);CHKERRQ(ierr);
60   ierr = PetscFree((*c)->queue);CHKERRQ(ierr);
61   ierr = PetscFree((*c)->queueLocal);CHKERRQ(ierr);
62   ierr = PetscFree((*c)->queueRemote);CHKERRQ(ierr);
63   ierr = PetscFree((*c)->neighbors);CHKERRQ(ierr);
64   ierr = PetscFree((*c)->needCount);CHKERRQ(ierr);
65   ierr = PetscFree((*c)->localOffsets);CHKERRQ(ierr);
66   ierr = PetscFree((*c)->fillCount);CHKERRQ(ierr);
67   ierr = PetscFree((*c)->remoteOffsets);CHKERRQ(ierr);
68   ierr = PetscFree((*c)->request);CHKERRQ(ierr);
69   ierr = PetscFree((*c)->status);CHKERRQ(ierr);
70   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "CharacteristicCreate"
76 PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
77 {
78   Characteristic newC;
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   PetscValidPointer(c, 2);
83   *c = PETSC_NULL;
84 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
85   ierr = CharacteristicInitializePackage(PETSC_NULL);CHKERRQ(ierr);
86 #endif
87 
88   ierr = PetscHeaderCreate(newC, _p_Characteristic, struct _CharacteristicOps, CHARACTERISTIC_CLASSID, -1, "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      = PETSC_NULL;
95   newC->velocity        = PETSC_NULL;
96   newC->velocityOld     = PETSC_NULL;
97   newC->numVelocityComp = 0;
98   newC->velocityComp    = PETSC_NULL;
99   newC->velocityInterp  = PETSC_NULL;
100   newC->velocityInterpLocal = PETSC_NULL;
101   newC->velocityCtx     = PETSC_NULL;
102   newC->fieldDA         = PETSC_NULL;
103   newC->field           = PETSC_NULL;
104   newC->numFieldComp    = 0;
105   newC->fieldComp       = PETSC_NULL;
106   newC->fieldInterp     = PETSC_NULL;
107   newC->fieldInterpLocal    = PETSC_NULL;
108   newC->fieldCtx        = PETSC_NULL;
109   newC->itemType        = PETSC_NULL;
110   newC->queue           = PETSC_NULL;
111   newC->queueSize       = 0;
112   newC->queueMax        = 0;
113   newC->queueLocal      = PETSC_NULL;
114   newC->queueLocalSize  = 0;
115   newC->queueLocalMax   = 0;
116   newC->queueRemote     = PETSC_NULL;
117   newC->queueRemoteSize = 0;
118   newC->queueRemoteMax  = 0;
119   newC->numNeighbors    = 0;
120   newC->neighbors       = PETSC_NULL;
121   newC->needCount       = PETSC_NULL;
122   newC->localOffsets    = PETSC_NULL;
123   newC->fillCount       = PETSC_NULL;
124   newC->remoteOffsets   = PETSC_NULL;
125   newC->request         = PETSC_NULL;
126   newC->status          = PETSC_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 = PETSC_NULL;
183     c->data = 0;
184   }
185 
186   ierr =  PetscFunctionListFind(((PetscObject)c)->comm,CharacteristicList, type,PETSC_TRUE, (void (**)(void)) &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 - See CharacteristicRegisterDynamic()
237 
238   Level: advanced
239 @*/
240 PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic))
241 {
242   PetscErrorCode ierr;
243   char           fullname[PETSC_MAX_PATH_LEN];
244 
245   PetscFunctionBegin;
246   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
247   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "CharacteristicSetVelocityInterpolation"
253 PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
254 {
255   PetscFunctionBegin;
256   c->velocityDA      = da;
257   c->velocity        = v;
258   c->velocityOld     = vOld;
259   c->numVelocityComp = numComponents;
260   c->velocityComp    = components;
261   c->velocityInterp  = interp;
262   c->velocityCtx     = ctx;
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal"
268 PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
269 {
270   PetscFunctionBegin;
271   c->velocityDA          = da;
272   c->velocity            = v;
273   c->velocityOld         = vOld;
274   c->numVelocityComp     = numComponents;
275   c->velocityComp        = components;
276   c->velocityInterpLocal = interp;
277   c->velocityCtx         = ctx;
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "CharacteristicSetFieldInterpolation"
283 PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
284 {
285   PetscFunctionBegin;
286 #if 0
287   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.");
288 #endif
289   c->fieldDA      = da;
290   c->field        = v;
291   c->numFieldComp = numComponents;
292   c->fieldComp    = components;
293   c->fieldInterp  = interp;
294   c->fieldCtx     = ctx;
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal"
300 PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal [], PetscInt, PetscInt [], PetscScalar [], void *), void *ctx)
301 {
302   PetscFunctionBegin;
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   PetscFunctionReturn(0);
313 }
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "CharacteristicSolve"
317 PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
318 {
319   CharacteristicPointDA2D Qi;
320   DM                      da = c->velocityDA;
321   Vec                     velocityLocal, velocityLocalOld;
322   Vec                     fieldLocal;
323   DMDALocalInfo           info;
324   PetscScalar             **solArray;
325   void                    *velocityArray;
326   void                    *velocityArrayOld;
327   void                    *fieldArray;
328   PassiveScalar           *interpIndices;
329   PassiveScalar           *velocityValues, *velocityValuesOld;
330   PassiveScalar           *fieldValues;
331   PetscMPIInt             rank;
332   PetscInt                dim;
333   PetscMPIInt             neighbors[9];
334   PetscInt                dof;
335   PetscInt                gx, gy;
336   PetscInt                n, is, ie, js, je, comp;
337   PetscErrorCode          ierr;
338 
339   PetscFunctionBegin;
340   c->queueSize = 0;
341   ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr);
342   ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr);
343   ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr);
344   ierr = CharacteristicSetUp(c);CHKERRQ(ierr);
345   /* global and local grid info */
346   ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr);
347   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
348   is   = info.xs;          ie   = info.xs+info.xm;
349   js   = info.ys;          je   = info.ys+info.ym;
350   /* Allocation */
351   ierr = PetscMalloc(dim*sizeof(PetscScalar),                &interpIndices);CHKERRQ(ierr);
352   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr);
353   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr);
354   ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar),    &fieldValues);CHKERRQ(ierr);
355   ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
356 
357   /* -----------------------------------------------------------------------
358      PART 1, AT t-dt/2
359      -----------------------------------------------------------------------*/
360   ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
361   /* GET POSITION AT HALF TIME IN THE PAST */
362   if (c->velocityInterpLocal) {
363     ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
364     ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
365     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
366     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
367     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
368     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
369     ierr = DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
370     ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
371   }
372   ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr);
373   for (Qi.j = js; Qi.j < je; Qi.j++) {
374     for (Qi.i = is; Qi.i < ie; Qi.i++) {
375       interpIndices[0] = Qi.i;
376       interpIndices[1] = Qi.j;
377       if (c->velocityInterpLocal) {
378         c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
379       } else {
380         c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
381       }
382       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
383       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
384 
385       /* Determine whether the position at t - dt/2 is local */
386       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
387 
388       /* Check for Periodic boundaries and move all periodic points back onto the domain */
389       ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
390       ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr);
391     }
392   }
393   ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
394 
395   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
396   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
397   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
398 
399   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
400   /* Calculate velocity at t_n+1/2 (local values) */
401   ierr = PetscInfo(PETSC_NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr);
402   for (n = 0; n < c->queueSize; n++) {
403     Qi = c->queue[n];
404     if (c->neighbors[Qi.proc] == rank) {
405       interpIndices[0] = Qi.x;
406       interpIndices[1] = Qi.y;
407       if (c->velocityInterpLocal) {
408         c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
409         c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
410       } else {
411         c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
412         c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
413       }
414       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
415       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
416     }
417     c->queue[n] = Qi;
418   }
419   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
420 
421   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
422   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
423   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
424 
425 
426   /* Calculate velocity at t_n+1/2 (fill remote requests) */
427   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
428   ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr);
429   for (n = 0; n < c->queueRemoteSize; n++) {
430     Qi = c->queueRemote[n];
431     interpIndices[0] = Qi.x;
432     interpIndices[1] = Qi.y;
433     if (c->velocityInterpLocal) {
434       c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
435       c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
436     } else {
437       c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
438       c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
439     }
440     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
441     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
442     c->queueRemote[n] = Qi;
443   }
444   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
445   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
446   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
447   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
448   if (c->velocityInterpLocal) {
449     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
450     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
451     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
452     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
453   }
454   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
455 
456   /* -----------------------------------------------------------------------
457      PART 2, AT t-dt
458      -----------------------------------------------------------------------*/
459 
460   /* GET POSITION AT t_n (local values) */
461   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
462   ierr = PetscInfo(PETSC_NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr);
463   for (n = 0; n < c->queueSize; n++) {
464     Qi = c->queue[n];
465     Qi.x = Qi.i - Qi.x*dt;
466     Qi.y = Qi.j - Qi.y*dt;
467 
468     /* Determine whether the position at t-dt is local */
469     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
470 
471     /* Check for Periodic boundaries and move all periodic points back onto the domain */
472     ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
473 
474     c->queue[n] = Qi;
475   }
476   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
477 
478   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
479   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
480   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
481 
482   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
483   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
484   if (c->fieldInterpLocal) {
485     ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
486     ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
487     ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
488     ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
489   }
490   ierr = PetscInfo(PETSC_NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr);
491   for (n = 0; n < c->queueSize; n++) {
492     if (c->neighbors[c->queue[n].proc] == rank) {
493       interpIndices[0] = c->queue[n].x;
494       interpIndices[1] = c->queue[n].y;
495       if (c->fieldInterpLocal) {
496         c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
497       } else {
498         c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
499       }
500       for (comp = 0; comp < c->numFieldComp; comp++) {
501         c->queue[n].field[comp] = fieldValues[comp];
502       }
503     }
504   }
505   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
506 
507   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
508   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
509   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
510 
511   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
512   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
513   ierr = PetscInfo1(PETSC_NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr);
514   for (n = 0; n < c->queueRemoteSize; n++) {
515     interpIndices[0] = c->queueRemote[n].x;
516     interpIndices[1] = c->queueRemote[n].y;
517 
518     /* for debugging purposes */
519     if (1) { /* hacked bounds test...let's do better */
520       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
521 
522       if ((im < (PetscScalar) is - 1.) || (im > (PetscScalar) ie) || (jm < (PetscScalar)  js - 1.) || (jm > (PetscScalar) je)) {
523         SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", im, jm);
524       }
525     }
526 
527     if (c->fieldInterpLocal) {
528       c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
529     } else {
530       c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
531     }
532     for (comp = 0; comp < c->numFieldComp; comp++) {
533       c->queueRemote[n].field[comp] = fieldValues[comp];
534     }
535   }
536   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
537 
538   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
539   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
540   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
541   if (c->fieldInterpLocal) {
542     ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
543     ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
544   }
545   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
546 
547   /* Return field of characteristics at t_n-1 */
548   ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
549   ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
550   ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
551   for (n = 0; n < c->queueSize; n++) {
552     Qi = c->queue[n];
553     for (comp = 0; comp < c->numFieldComp; comp++) {
554       solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
555     }
556   }
557   ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
558   ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
559   ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
560 
561   /* Cleanup */
562   ierr = PetscFree(interpIndices);CHKERRQ(ierr);
563   ierr = PetscFree(velocityValues);CHKERRQ(ierr);
564   ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr);
565   ierr = PetscFree(fieldValues);CHKERRQ(ierr);
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "CharacteristicSetNeighbors"
571 PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
572 {
573   PetscErrorCode ierr;
574 
575   PetscFunctionBegin;
576   c->numNeighbors = numNeighbors;
577   ierr = PetscFree(c->neighbors);CHKERRQ(ierr);
578   ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr);
579   ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr);
580   PetscFunctionReturn(0);
581 }
582 
583 #undef __FUNCT__
584 #define __FUNCT__ "CharacteristicAddPoint"
585 PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
586 {
587   PetscFunctionBegin;
588   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
589   c->queue[c->queueSize++] = *point;
590   PetscFunctionReturn(0);
591 }
592 
593 #undef __FUNCT__
594 #define __FUNCT__ "CharacteristicSendCoordinatesBegin"
595 int CharacteristicSendCoordinatesBegin(Characteristic c)
596 {
597   PetscMPIInt    rank, tag = 121;
598   PetscInt       i, n;
599   PetscErrorCode ierr;
600 
601   PetscFunctionBegin;
602   ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr);
603   ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr);
604   ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr);
605   for (i = 0;  i < c->queueSize; i++) {
606     c->needCount[c->queue[i].proc]++;
607   }
608   c->fillCount[0] = 0;
609   for (n = 1; n < c->numNeighbors; n++) {
610     ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr);
611   }
612   for (n = 1; n < c->numNeighbors; n++) {
613     ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr);
614   }
615   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
616   /* Initialize the remote queue */
617   c->queueLocalMax  = c->localOffsets[0]  = 0;
618   c->queueRemoteMax = c->remoteOffsets[0] = 0;
619   for (n = 1; n < c->numNeighbors; n++) {
620     c->remoteOffsets[n] = c->queueRemoteMax;
621     c->queueRemoteMax  += c->fillCount[n];
622     c->localOffsets[n]  = c->queueLocalMax;
623     c->queueLocalMax   += c->needCount[n];
624   }
625   /* HACK BEGIN */
626   for (n = 1; n < c->numNeighbors; n++) {
627     c->localOffsets[n] += c->needCount[0];
628   }
629   c->needCount[0] = 0;
630   /* HACK END */
631   if (c->queueRemoteMax) {
632     ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr);
633   } else {
634     c->queueRemote = PETSC_NULL;
635   }
636   c->queueRemoteSize = c->queueRemoteMax;
637 
638   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
639   for (n = 1; n < c->numNeighbors; n++) {
640     ierr = PetscInfo2(PETSC_NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr);
641     ierr = MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr);
642   }
643   for (n = 1; n < c->numNeighbors; n++) {
644     ierr = PetscInfo2(PETSC_NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr);
645     ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr);
646   }
647   PetscFunctionReturn(0);
648 }
649 
650 #undef __FUNCT__
651 #define __FUNCT__ "CharacteristicSendCoordinatesEnd"
652 PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
653 {
654 #if 0
655   PetscMPIInt rank;
656   PetscInt n;
657 #endif
658   PetscErrorCode ierr;
659 
660   PetscFunctionBegin;
661   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
662 #if 0
663   ierr = MPI_Comm_rank(((PetscObject)c)->comm, &rank);CHKERRQ(ierr);
664   for (n = 0; n < c->queueRemoteSize; n++) {
665     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);
666   }
667 #endif
668   PetscFunctionReturn(0);
669 }
670 
671 #undef __FUNCT__
672 #define __FUNCT__ "CharacteristicGetValuesBegin"
673 PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
674 {
675   PetscMPIInt    tag = 121;
676   PetscInt       n;
677   PetscErrorCode ierr;
678 
679   PetscFunctionBegin;
680   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
681   for (n = 1; n < c->numNeighbors; n++) {
682     ierr = MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm, &(c->request[n-1]));CHKERRQ(ierr);
683   }
684   for (n = 1; n < c->numNeighbors; n++) {
685     ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, ((PetscObject)c)->comm);CHKERRQ(ierr);
686   }
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "CharacteristicGetValuesEnd"
692 PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
693 {
694   PetscErrorCode ierr;
695 
696   PetscFunctionBegin;
697   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
698   /* Free queue of requests from other procs */
699   ierr = PetscFree(c->queueRemote);CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 
703 /*---------------------------------------------------------------------*/
704 #undef __FUNCT__
705 #define __FUNCT__ "CharacteristicHeapSort"
706 /*
707   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
708 */
709 PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
710 /*---------------------------------------------------------------------*/
711 {
712   PetscErrorCode          ierr;
713   CharacteristicPointDA2D temp;
714   PetscInt                n;
715 
716   PetscFunctionBegin;
717   if (0) { /* Check the order of the queue before sorting */
718     PetscInfo(PETSC_NULL, "Before Heap sort\n");
719     for (n=0;  n<size; n++) {
720       ierr = PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
721     }
722   }
723 
724   /* SORTING PHASE */
725   for (n = (size / 2)-1; n >= 0; n--) {
726     ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/
727   }
728   for (n = size-1; n >= 1; n--) {
729     temp = queue[0];
730     queue[0] = queue[n];
731     queue[n] = temp;
732     ierr = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr);
733   }
734   if (0) { /* Check the order of the queue after sorting */
735     ierr = PetscInfo(PETSC_NULL, "Avter  Heap sort\n");CHKERRQ(ierr);
736     for (n=0;  n<size; n++) {
737       ierr = PetscInfo2(PETSC_NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
738     }
739   }
740   PetscFunctionReturn(0);
741 }
742 
743 /*---------------------------------------------------------------------*/
744 #undef __FUNCT__
745 #define __FUNCT__ "CharacteristicSiftDown"
746 /*
747   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
748 */
749 PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
750 /*---------------------------------------------------------------------*/
751 {
752   PetscBool                done = PETSC_FALSE;
753   PetscInt                 maxChild;
754   CharacteristicPointDA2D  temp;
755 
756   PetscFunctionBegin;
757   while ((root*2 <= bottom) && (!done)) {
758     if (root*2 == bottom)  maxChild = root * 2;
759     else if (queue[root*2].proc > queue[root*2+1].proc)  maxChild = root * 2;
760     else  maxChild = root * 2 + 1;
761 
762     if (queue[root].proc < queue[maxChild].proc) {
763       temp = queue[root];
764       queue[root] = queue[maxChild];
765       queue[maxChild] = temp;
766       root = maxChild;
767     } else done = PETSC_TRUE;
768   }
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "DMDAGetNeighborsRank"
774 /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
775 PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
776 {
777   DMDABoundaryType bx, by;
778   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
779   MPI_Comm         comm;
780   PetscMPIInt      rank;
781   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
782   PetscErrorCode   ierr;
783 
784   PetscFunctionBegin;
785   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
786   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
787   ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr);
788 
789   if (bx == DMDA_BOUNDARY_PERIODIC) {
790     IPeriodic = PETSC_TRUE;
791   }
792   if (by == DMDA_BOUNDARY_PERIODIC) {
793     JPeriodic = PETSC_TRUE;
794   }
795 
796   neighbors[0] = rank;
797   rank = 0;
798   ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr);
799   for (pj=0;pj<PJ;pj++) {
800     ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr);
801     for (pi=0;pi<PI;pi++) {
802       procs[pj][pi] = rank;
803       rank++;
804     }
805   }
806 
807   pi = neighbors[0] % PI;
808   pj = neighbors[0] / PI;
809   pim = pi-1;  if (pim<0) pim=PI-1;
810   pip = (pi+1)%PI;
811   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
812   pjp = (pj+1)%PJ;
813 
814   neighbors[1] = procs[pj] [pim];
815   neighbors[2] = procs[pjp][pim];
816   neighbors[3] = procs[pjp][pi];
817   neighbors[4] = procs[pjp][pip];
818   neighbors[5] = procs[pj] [pip];
819   neighbors[6] = procs[pjm][pip];
820   neighbors[7] = procs[pjm][pi];
821   neighbors[8] = procs[pjm][pim];
822 
823   if (!IPeriodic) {
824     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
825     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
826   }
827 
828   if (!JPeriodic) {
829     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
830     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
831   }
832 
833   for (pj = 0; pj < PJ; pj++) {
834     ierr = PetscFree(procs[pj]);CHKERRQ(ierr);
835   }
836   ierr = PetscFree(procs);CHKERRQ(ierr);
837   PetscFunctionReturn(0);
838 }
839 
840 #undef __FUNCT__
841 #define __FUNCT__ "DMDAGetNeighborRelative"
842 /*
843   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
844     2 | 3 | 4
845     __|___|__
846     1 | 0 | 5
847     __|___|__
848     8 | 7 | 6
849       |   |
850 */
851 PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr)
852 {
853   DMDALocalInfo  info;
854   PassiveReal    is,ie,js,je;
855   PetscErrorCode ierr;
856 
857   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
858   is = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5;
859   js = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5;
860 
861   if (ir >= is && ir <= ie) { /* center column */
862     if (jr >= js && jr <= je) {
863       return 0;
864     } else if (jr < js) {
865       return 7;
866     } else {
867       return 3;
868     }
869   } else if (ir < is) {     /* left column */
870     if (jr >= js && jr <= je) {
871       return 1;
872     } else if (jr < js) {
873       return 8;
874     } else {
875       return 2;
876     }
877   } else {                  /* right column */
878     if (jr >= js && jr <= je) {
879       return 5;
880     } else if (jr < js) {
881       return 6;
882     } else {
883       return 4;
884     }
885   }
886 }
887