xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 1c84c290429673533d0af1ef36b216eb00726fbb)
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 -  Adds a solver to the method of characteristics package.
239 
240    Not Collective
241 
242    Input Parameters:
243 +  name_solver - name of a new user-defined solver
244 .  path - path (either absolute or relative) the library containing this solver
245 .  name_create - name of routine to create method context
246 -  routine_create - routine to create method context
247 
248   Sample usage:
249 .vb
250     CharacteristicRegister("my_char","MyCharCreate", MyCharCreate);
251 .ve
252 
253   Then, your Characteristic type can be chosen with the procedural interface via
254 .vb
255     CharacteristicCreate(MPI_Comm, Characteristic* &char);
256     CharacteristicSetType(char,"my_char");
257 .ve
258    or at runtime via the option
259 .vb
260     -characteristic_type my_char
261 .ve
262 
263    Notes:
264    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
265 
266 .keywords: Characteristic, register
267 
268 .seealso: CharacteristicRegisterAll(), CharacteristicRegisterDestroy()
269 
270   Level: advanced
271 @*/
272 PetscErrorCode CharacteristicRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Characteristic))
273 {
274   PetscErrorCode ierr;
275   char           fullname[PETSC_MAX_PATH_LEN];
276 
277   PetscFunctionBegin;
278   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
279   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&CharacteristicList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "CharacteristicSetVelocityInterpolation"
285 PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
286 {
287   PetscFunctionBegin;
288   c->velocityDA      = da;
289   c->velocity        = v;
290   c->velocityOld     = vOld;
291   c->numVelocityComp = numComponents;
292   c->velocityComp    = components;
293   c->velocityInterp  = interp;
294   c->velocityCtx     = ctx;
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "CharacteristicSetVelocityInterpolationLocal"
300 PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal [], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
301 {
302   PetscFunctionBegin;
303   c->velocityDA          = da;
304   c->velocity            = v;
305   c->velocityOld         = vOld;
306   c->numVelocityComp     = numComponents;
307   c->velocityComp        = components;
308   c->velocityInterpLocal = interp;
309   c->velocityCtx         = ctx;
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "CharacteristicSetFieldInterpolation"
315 PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
316 {
317   PetscFunctionBegin;
318 #if 0
319   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.");
320 #endif
321   c->fieldDA      = da;
322   c->field        = v;
323   c->numFieldComp = numComponents;
324   c->fieldComp    = components;
325   c->fieldInterp  = interp;
326   c->fieldCtx     = ctx;
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "CharacteristicSetFieldInterpolationLocal"
332 PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
333 {
334   PetscFunctionBegin;
335 #if 0
336   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.");
337 #endif
338   c->fieldDA          = da;
339   c->field            = v;
340   c->numFieldComp     = numComponents;
341   c->fieldComp        = components;
342   c->fieldInterpLocal = interp;
343   c->fieldCtx         = ctx;
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "CharacteristicSolve"
349 PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
350 {
351   CharacteristicPointDA2D Qi;
352   DM                      da = c->velocityDA;
353   Vec                     velocityLocal, velocityLocalOld;
354   Vec                     fieldLocal;
355   DMDALocalInfo           info;
356   PetscScalar             **solArray;
357   void                    *velocityArray;
358   void                    *velocityArrayOld;
359   void                    *fieldArray;
360   PassiveScalar           *interpIndices;
361   PassiveScalar           *velocityValues, *velocityValuesOld;
362   PassiveScalar           *fieldValues;
363   PetscMPIInt             rank;
364   PetscInt                dim;
365   PetscMPIInt             neighbors[9];
366   PetscInt                dof;
367   PetscInt                gx, gy;
368   PetscInt                n, is, ie, js, je, comp;
369   PetscErrorCode          ierr;
370 
371   PetscFunctionBegin;
372   c->queueSize = 0;
373   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
374   ierr = DMDAGetNeighborsRank(da, neighbors);CHKERRQ(ierr);
375   ierr = CharacteristicSetNeighbors(c, 9, neighbors);CHKERRQ(ierr);
376   ierr = CharacteristicSetUp(c);CHKERRQ(ierr);
377   /* global and local grid info */
378   ierr = DMDAGetInfo(da, &dim, &gx, &gy, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0);CHKERRQ(ierr);
379   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
380   is   = info.xs;          ie   = info.xs+info.xm;
381   js   = info.ys;          je   = info.ys+info.ym;
382   /* Allocation */
383   ierr = PetscMalloc(dim*sizeof(PetscScalar),                &interpIndices);CHKERRQ(ierr);
384   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValues);CHKERRQ(ierr);
385   ierr = PetscMalloc(c->numVelocityComp*sizeof(PetscScalar), &velocityValuesOld);CHKERRQ(ierr);
386   ierr = PetscMalloc(c->numFieldComp*sizeof(PetscScalar),    &fieldValues);CHKERRQ(ierr);
387   ierr = PetscLogEventBegin(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
388 
389   /* -----------------------------------------------------------------------
390      PART 1, AT t-dt/2
391      -----------------------------------------------------------------------*/
392   ierr = PetscLogEventBegin(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
393   /* GET POSITION AT HALF TIME IN THE PAST */
394   if (c->velocityInterpLocal) {
395     ierr = DMGetLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
396     ierr = DMGetLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
397     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
398     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal);CHKERRQ(ierr);
399     ierr = DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
400     ierr = DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld);CHKERRQ(ierr);
401     ierr = DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
402     ierr = DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
403   }
404   ierr = PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n");CHKERRQ(ierr);
405   for (Qi.j = js; Qi.j < je; Qi.j++) {
406     for (Qi.i = is; Qi.i < ie; Qi.i++) {
407       interpIndices[0] = Qi.i;
408       interpIndices[1] = Qi.j;
409       if (c->velocityInterpLocal) c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
410       else c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx);
411       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
412       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
413 
414       /* Determine whether the position at t - dt/2 is local */
415       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
416 
417       /* Check for Periodic boundaries and move all periodic points back onto the domain */
418       ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
419       ierr = CharacteristicAddPoint(c, &Qi);CHKERRQ(ierr);
420     }
421   }
422   ierr = PetscLogEventEnd(CHARACTERISTIC_QueueSetup,0,0,0,0);CHKERRQ(ierr);
423 
424   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
425   ierr = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
426   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
427 
428   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
429   /* Calculate velocity at t_n+1/2 (local values) */
430   ierr = PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n");CHKERRQ(ierr);
431   for (n = 0; n < c->queueSize; n++) {
432     Qi = c->queue[n];
433     if (c->neighbors[Qi.proc] == rank) {
434       interpIndices[0] = Qi.x;
435       interpIndices[1] = Qi.y;
436       if (c->velocityInterpLocal) {
437         c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
438         c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
439       } else {
440         c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
441         c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
442       }
443       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
444       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
445     }
446     c->queue[n] = Qi;
447   }
448   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,0,0,0,0);CHKERRQ(ierr);
449 
450   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
451   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
452   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
453 
454 
455   /* Calculate velocity at t_n+1/2 (fill remote requests) */
456   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
457   ierr = PetscInfo1(NULL, "Calculating %d remote velocities at t_{n - 1/2}\n", c->queueRemoteSize);CHKERRQ(ierr);
458   for (n = 0; n < c->queueRemoteSize; n++) {
459     Qi = c->queueRemote[n];
460     interpIndices[0] = Qi.x;
461     interpIndices[1] = Qi.y;
462     if (c->velocityInterpLocal) {
463       c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
464       c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
465     } else {
466       c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx);
467       c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx);
468     }
469     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
470     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
471     c->queueRemote[n] = Qi;
472   }
473   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,0,0,0,0);CHKERRQ(ierr);
474   ierr = PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
475   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
476   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
477   if (c->velocityInterpLocal) {
478     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray);CHKERRQ(ierr);
479     ierr = DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld);CHKERRQ(ierr);
480     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocal);CHKERRQ(ierr);
481     ierr = DMRestoreLocalVector(c->velocityDA, &velocityLocalOld);CHKERRQ(ierr);
482   }
483   ierr = PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,0,0,0,0);CHKERRQ(ierr);
484 
485   /* -----------------------------------------------------------------------
486      PART 2, AT t-dt
487      -----------------------------------------------------------------------*/
488 
489   /* GET POSITION AT t_n (local values) */
490   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
491   ierr = PetscInfo(NULL, "Calculating position at t_{n}\n");CHKERRQ(ierr);
492   for (n = 0; n < c->queueSize; n++) {
493     Qi   = c->queue[n];
494     Qi.x = Qi.i - Qi.x*dt;
495     Qi.y = Qi.j - Qi.y*dt;
496 
497     /* Determine whether the position at t-dt is local */
498     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
499 
500     /* Check for Periodic boundaries and move all periodic points back onto the domain */
501     ierr = DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y));CHKERRQ(ierr);
502 
503     c->queue[n] = Qi;
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 = CharacteristicSendCoordinatesBegin(c);CHKERRQ(ierr);
509   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
510 
511   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
512   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
513   if (c->fieldInterpLocal) {
514     ierr = DMGetLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
515     ierr = DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
516     ierr = DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal);CHKERRQ(ierr);
517     ierr = DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
518   }
519   ierr = PetscInfo(NULL, "Calculating local field at t_{n}\n");CHKERRQ(ierr);
520   for (n = 0; n < c->queueSize; n++) {
521     if (c->neighbors[c->queue[n].proc] == rank) {
522       interpIndices[0] = c->queue[n].x;
523       interpIndices[1] = c->queue[n].y;
524       if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
525       else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
526       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
527     }
528   }
529   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,0,0,0,0);CHKERRQ(ierr);
530 
531   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
532   ierr = CharacteristicSendCoordinatesEnd(c);CHKERRQ(ierr);
533   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
534 
535   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
536   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
537   ierr = PetscInfo1(NULL, "Calculating %d remote field points at t_{n}\n", c->queueRemoteSize);CHKERRQ(ierr);
538   for (n = 0; n < c->queueRemoteSize; n++) {
539     interpIndices[0] = c->queueRemote[n].x;
540     interpIndices[1] = c->queueRemote[n].y;
541 
542     /* for debugging purposes */
543     if (1) { /* hacked bounds test...let's do better */
544       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
545 
546       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);
547     }
548 
549     if (c->fieldInterpLocal) c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
550     else c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx);
551     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
552   }
553   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,0,0,0,0);CHKERRQ(ierr);
554 
555   ierr = PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
556   ierr = CharacteristicGetValuesBegin(c);CHKERRQ(ierr);
557   ierr = CharacteristicGetValuesEnd(c);CHKERRQ(ierr);
558   if (c->fieldInterpLocal) {
559     ierr = DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray);CHKERRQ(ierr);
560     ierr = DMRestoreLocalVector(c->fieldDA, &fieldLocal);CHKERRQ(ierr);
561   }
562   ierr = PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,0,0,0,0);CHKERRQ(ierr);
563 
564   /* Return field of characteristics at t_n-1 */
565   ierr = PetscLogEventBegin(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
566   ierr = DMDAGetInfo(c->fieldDA,0,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
567   ierr = DMDAVecGetArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
568   for (n = 0; n < c->queueSize; n++) {
569     Qi = c->queue[n];
570     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
571   }
572   ierr = DMDAVecRestoreArray(c->fieldDA, solution, &solArray);CHKERRQ(ierr);
573   ierr = PetscLogEventEnd(CHARACTERISTIC_DAUpdate,0,0,0,0);CHKERRQ(ierr);
574   ierr = PetscLogEventEnd(CHARACTERISTIC_Solve,0,0,0,0);CHKERRQ(ierr);
575 
576   /* Cleanup */
577   ierr = PetscFree(interpIndices);CHKERRQ(ierr);
578   ierr = PetscFree(velocityValues);CHKERRQ(ierr);
579   ierr = PetscFree(velocityValuesOld);CHKERRQ(ierr);
580   ierr = PetscFree(fieldValues);CHKERRQ(ierr);
581   PetscFunctionReturn(0);
582 }
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "CharacteristicSetNeighbors"
586 PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
587 {
588   PetscErrorCode ierr;
589 
590   PetscFunctionBegin;
591   c->numNeighbors = numNeighbors;
592   ierr = PetscFree(c->neighbors);CHKERRQ(ierr);
593   ierr = PetscMalloc(numNeighbors * sizeof(PetscMPIInt), &c->neighbors);CHKERRQ(ierr);
594   ierr = PetscMemcpy(c->neighbors, neighbors, numNeighbors * sizeof(PetscMPIInt));CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "CharacteristicAddPoint"
600 PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
601 {
602   PetscFunctionBegin;
603   if (c->queueSize >= c->queueMax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeeded maximum queue size %d", c->queueMax);
604   c->queue[c->queueSize++] = *point;
605   PetscFunctionReturn(0);
606 }
607 
608 #undef __FUNCT__
609 #define __FUNCT__ "CharacteristicSendCoordinatesBegin"
610 int CharacteristicSendCoordinatesBegin(Characteristic c)
611 {
612   PetscMPIInt    rank, tag = 121;
613   PetscInt       i, n;
614   PetscErrorCode ierr;
615 
616   PetscFunctionBegin;
617   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
618   ierr = CharacteristicHeapSort(c, c->queue, c->queueSize);CHKERRQ(ierr);
619   ierr = PetscMemzero(c->needCount, c->numNeighbors * sizeof(PetscInt));CHKERRQ(ierr);
620   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
621   c->fillCount[0] = 0;
622   for (n = 1; n < c->numNeighbors; n++) {
623     ierr = MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1]));CHKERRQ(ierr);
624   }
625   for (n = 1; n < c->numNeighbors; n++) {
626     ierr = MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
627   }
628   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
629   /* Initialize the remote queue */
630   c->queueLocalMax  = c->localOffsets[0]  = 0;
631   c->queueRemoteMax = c->remoteOffsets[0] = 0;
632   for (n = 1; n < c->numNeighbors; n++) {
633     c->remoteOffsets[n] = c->queueRemoteMax;
634     c->queueRemoteMax  += c->fillCount[n];
635     c->localOffsets[n]  = c->queueLocalMax;
636     c->queueLocalMax   += c->needCount[n];
637   }
638   /* HACK BEGIN */
639   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
640   c->needCount[0] = 0;
641   /* HACK END */
642   if (c->queueRemoteMax) {
643     ierr = PetscMalloc(sizeof(CharacteristicPointDA2D) * c->queueRemoteMax, &c->queueRemote);CHKERRQ(ierr);
644   } else c->queueRemote = NULL;
645   c->queueRemoteSize = c->queueRemoteMax;
646 
647   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
648   for (n = 1; n < c->numNeighbors; n++) {
649     ierr = PetscInfo2(NULL, "Receiving %d requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]);CHKERRQ(ierr);
650     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);
651   }
652   for (n = 1; n < c->numNeighbors; n++) {
653     ierr = PetscInfo2(NULL, "Sending %d requests for values from proc %d\n", c->needCount[n], c->neighbors[n]);CHKERRQ(ierr);
654     ierr = MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
655   }
656   PetscFunctionReturn(0);
657 }
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "CharacteristicSendCoordinatesEnd"
661 PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
662 {
663 #if 0
664   PetscMPIInt rank;
665   PetscInt    n;
666 #endif
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
671 #if 0
672   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank);CHKERRQ(ierr);
673   for (n = 0; n < c->queueRemoteSize; n++) {
674     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);
675   }
676 #endif
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "CharacteristicGetValuesBegin"
682 PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
683 {
684   PetscMPIInt    tag = 121;
685   PetscInt       n;
686   PetscErrorCode ierr;
687 
688   PetscFunctionBegin;
689   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
690   for (n = 1; n < c->numNeighbors; n++) {
691     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);
692   }
693   for (n = 1; n < c->numNeighbors; n++) {
694     ierr = MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c));CHKERRQ(ierr);
695   }
696   PetscFunctionReturn(0);
697 }
698 
699 #undef __FUNCT__
700 #define __FUNCT__ "CharacteristicGetValuesEnd"
701 PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   ierr = MPI_Waitall(c->numNeighbors-1, c->request, c->status);CHKERRQ(ierr);
707   /* Free queue of requests from other procs */
708   ierr = PetscFree(c->queueRemote);CHKERRQ(ierr);
709   PetscFunctionReturn(0);
710 }
711 
712 /*---------------------------------------------------------------------*/
713 #undef __FUNCT__
714 #define __FUNCT__ "CharacteristicHeapSort"
715 /*
716   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
717 */
718 PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
719 /*---------------------------------------------------------------------*/
720 {
721   PetscErrorCode          ierr;
722   CharacteristicPointDA2D temp;
723   PetscInt                n;
724 
725   PetscFunctionBegin;
726   if (0) { /* Check the order of the queue before sorting */
727     PetscInfo(NULL, "Before Heap sort\n");
728     for (n=0; n<size; n++) {
729       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
730     }
731   }
732 
733   /* SORTING PHASE */
734   for (n = (size / 2)-1; n >= 0; n--) {
735     ierr = CharacteristicSiftDown(c, queue, n, size-1);CHKERRQ(ierr); /* Rich had size-1 here, Matt had size*/
736   }
737   for (n = size-1; n >= 1; n--) {
738     temp     = queue[0];
739     queue[0] = queue[n];
740     queue[n] = temp;
741     ierr     = CharacteristicSiftDown(c, queue, 0, n-1);CHKERRQ(ierr);
742   }
743   if (0) { /* Check the order of the queue after sorting */
744     ierr = PetscInfo(NULL, "Avter  Heap sort\n");CHKERRQ(ierr);
745     for (n=0; n<size; n++) {
746       ierr = PetscInfo2(NULL,"%d %d\n",n,queue[n].proc);CHKERRQ(ierr);
747     }
748   }
749   PetscFunctionReturn(0);
750 }
751 
752 /*---------------------------------------------------------------------*/
753 #undef __FUNCT__
754 #define __FUNCT__ "CharacteristicSiftDown"
755 /*
756   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
757 */
758 PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
759 /*---------------------------------------------------------------------*/
760 {
761   PetscBool               done = PETSC_FALSE;
762   PetscInt                maxChild;
763   CharacteristicPointDA2D temp;
764 
765   PetscFunctionBegin;
766   while ((root*2 <= bottom) && (!done)) {
767     if (root*2 == bottom) maxChild = root * 2;
768     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
769     else maxChild = root * 2 + 1;
770 
771     if (queue[root].proc < queue[maxChild].proc) {
772       temp            = queue[root];
773       queue[root]     = queue[maxChild];
774       queue[maxChild] = temp;
775       root            = maxChild;
776     } else done = PETSC_TRUE;
777   }
778   PetscFunctionReturn(0);
779 }
780 
781 #undef __FUNCT__
782 #define __FUNCT__ "DMDAGetNeighborsRank"
783 /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
784 PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
785 {
786   DMDABoundaryType bx, by;
787   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
788   MPI_Comm         comm;
789   PetscMPIInt      rank;
790   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
791   PetscErrorCode   ierr;
792 
793   PetscFunctionBegin;
794   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
795   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
796   ierr = DMDAGetInfo(da, 0, 0, 0, 0, &PI,&PJ, 0, 0, 0, &bx, &by,0, 0);CHKERRQ(ierr);
797 
798   if (bx == DMDA_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
799   if (by == DMDA_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
800 
801   neighbors[0] = rank;
802   rank = 0;
803   ierr = PetscMalloc(sizeof(PetscInt*)*PJ,&procs);CHKERRQ(ierr);
804   for (pj=0; pj<PJ; pj++) {
805     ierr = PetscMalloc(sizeof(PetscInt)*PI,&(procs[pj]));CHKERRQ(ierr);
806     for (pi=0; pi<PI; pi++) {
807       procs[pj][pi] = rank;
808       rank++;
809     }
810   }
811 
812   pi  = neighbors[0] % PI;
813   pj  = neighbors[0] / PI;
814   pim = pi-1;  if (pim<0) pim=PI-1;
815   pip = (pi+1)%PI;
816   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
817   pjp = (pj+1)%PJ;
818 
819   neighbors[1] = procs[pj] [pim];
820   neighbors[2] = procs[pjp][pim];
821   neighbors[3] = procs[pjp][pi];
822   neighbors[4] = procs[pjp][pip];
823   neighbors[5] = procs[pj] [pip];
824   neighbors[6] = procs[pjm][pip];
825   neighbors[7] = procs[pjm][pi];
826   neighbors[8] = procs[pjm][pim];
827 
828   if (!IPeriodic) {
829     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
830     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
831   }
832 
833   if (!JPeriodic) {
834     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
835     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
836   }
837 
838   for (pj = 0; pj < PJ; pj++) {
839     ierr = PetscFree(procs[pj]);CHKERRQ(ierr);
840   }
841   ierr = PetscFree(procs);CHKERRQ(ierr);
842   PetscFunctionReturn(0);
843 }
844 
845 #undef __FUNCT__
846 #define __FUNCT__ "DMDAGetNeighborRelative"
847 /*
848   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
849     2 | 3 | 4
850     __|___|__
851     1 | 0 | 5
852     __|___|__
853     8 | 7 | 6
854       |   |
855 */
856 PetscInt DMDAGetNeighborRelative(DM da, PassiveReal ir, PassiveReal jr)
857 {
858   DMDALocalInfo  info;
859   PassiveReal    is,ie,js,je;
860   PetscErrorCode ierr;
861 
862   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
863   is   = (PassiveReal) info.xs - 0.5; ie = (PassiveReal) info.xs + info.xm - 0.5;
864   js   = (PassiveReal) info.ys - 0.5; je = (PassiveReal) info.ys + info.ym - 0.5;
865 
866   if (ir >= is && ir <= ie) { /* center column */
867     if (jr >= js && jr <= je) return 0;
868     else if (jr < js)         return 7;
869     else                      return 3;
870   } else if (ir < is) {     /* left column */
871     if (jr >= js && jr <= je) return 1;
872     else if (jr < js)         return 8;
873     else                      return 2;
874   } else {                  /* right column */
875     if (jr >= js && jr <= je) return 5;
876     else if (jr < js)         return 6;
877     else                      return 4;
878   }
879 }
880