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