xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1af33a6ddSJed Brown 
2af0996ceSBarry Smith #include <petsc/private/characteristicimpl.h> /*I "petsccharacteristic.h" I*/
31e25c274SJed Brown #include <petscdmda.h>
4665c2dedSJed Brown #include <petscviewer.h>
5af33a6ddSJed Brown 
6af33a6ddSJed Brown PetscClassId  CHARACTERISTIC_CLASSID;
7af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_SetUp, CHARACTERISTIC_Solve, CHARACTERISTIC_QueueSetup, CHARACTERISTIC_DAUpdate;
8af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_HalfTimeLocal, CHARACTERISTIC_HalfTimeRemote, CHARACTERISTIC_HalfTimeExchange;
9af33a6ddSJed Brown PetscLogEvent CHARACTERISTIC_FullTimeLocal, CHARACTERISTIC_FullTimeRemote, CHARACTERISTIC_FullTimeExchange;
10af33a6ddSJed Brown /*
11af33a6ddSJed Brown    Contains the list of registered characteristic routines
12af33a6ddSJed Brown */
130298fd71SBarry Smith PetscFunctionList CharacteristicList              = NULL;
14140e18c1SBarry Smith PetscBool         CharacteristicRegisterAllCalled = PETSC_FALSE;
15af33a6ddSJed Brown 
16af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt []);
171cc8b206SBarry Smith PetscInt       DMDAGetNeighborRelative(DM, PetscReal, PetscReal);
18af33a6ddSJed Brown PetscErrorCode DMDAMapToPeriodicDomain(DM, PetscScalar []);
19af33a6ddSJed Brown 
200b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
210b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
22af33a6ddSJed Brown 
23af33a6ddSJed Brown PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24af33a6ddSJed Brown {
25af33a6ddSJed Brown   PetscBool      iascii;
26af33a6ddSJed Brown 
27af33a6ddSJed Brown   PetscFunctionBegin;
28af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
29af33a6ddSJed Brown   if (!viewer) {
309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c),&viewer));
31af33a6ddSJed Brown   }
32af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
33af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
34af33a6ddSJed Brown 
359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
36af33a6ddSJed Brown   if (!iascii) {
37*dbbe0bcdSBarry Smith     PetscTryTypeMethod(c,view, viewer);
38af33a6ddSJed Brown   }
39af33a6ddSJed Brown   PetscFunctionReturn(0);
40af33a6ddSJed Brown }
41af33a6ddSJed Brown 
426bf464f9SBarry Smith PetscErrorCode CharacteristicDestroy(Characteristic *c)
43af33a6ddSJed Brown {
44af33a6ddSJed Brown   PetscFunctionBegin;
456bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
466bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
476bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
48af33a6ddSJed Brown 
496bf464f9SBarry Smith   if ((*c)->ops->destroy) {
509566063dSJacob Faibussowitsch     PetscCall((*(*c)->ops->destroy)((*c)));
51af33a6ddSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&(*c)->itemType));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queue));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueLocal));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueRemote));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->neighbors));
579566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->needCount));
589566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->localOffsets));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->fillCount));
609566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->remoteOffsets));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->request));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->status));
639566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
64af33a6ddSJed Brown   PetscFunctionReturn(0);
65af33a6ddSJed Brown }
66af33a6ddSJed Brown 
67af33a6ddSJed Brown PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
68af33a6ddSJed Brown {
69af33a6ddSJed Brown   Characteristic newC;
70af33a6ddSJed Brown 
71af33a6ddSJed Brown   PetscFunctionBegin;
72af33a6ddSJed Brown   PetscValidPointer(c, 2);
730298fd71SBarry Smith   *c = NULL;
749566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
75af33a6ddSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
77af33a6ddSJed Brown   *c   = newC;
78af33a6ddSJed Brown 
79af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
80af33a6ddSJed Brown   newC->numIds              = 0;
810298fd71SBarry Smith   newC->velocityDA          = NULL;
820298fd71SBarry Smith   newC->velocity            = NULL;
830298fd71SBarry Smith   newC->velocityOld         = NULL;
84af33a6ddSJed Brown   newC->numVelocityComp     = 0;
850298fd71SBarry Smith   newC->velocityComp        = NULL;
860298fd71SBarry Smith   newC->velocityInterp      = NULL;
870298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
880298fd71SBarry Smith   newC->velocityCtx         = NULL;
890298fd71SBarry Smith   newC->fieldDA             = NULL;
900298fd71SBarry Smith   newC->field               = NULL;
91af33a6ddSJed Brown   newC->numFieldComp        = 0;
920298fd71SBarry Smith   newC->fieldComp           = NULL;
930298fd71SBarry Smith   newC->fieldInterp         = NULL;
940298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
950298fd71SBarry Smith   newC->fieldCtx            = NULL;
960298fd71SBarry Smith   newC->itemType            = 0;
970298fd71SBarry Smith   newC->queue               = NULL;
98af33a6ddSJed Brown   newC->queueSize           = 0;
99af33a6ddSJed Brown   newC->queueMax            = 0;
1000298fd71SBarry Smith   newC->queueLocal          = NULL;
101af33a6ddSJed Brown   newC->queueLocalSize      = 0;
102af33a6ddSJed Brown   newC->queueLocalMax       = 0;
1030298fd71SBarry Smith   newC->queueRemote         = NULL;
104af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
105af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
106af33a6ddSJed Brown   newC->numNeighbors        = 0;
1070298fd71SBarry Smith   newC->neighbors           = NULL;
1080298fd71SBarry Smith   newC->needCount           = NULL;
1090298fd71SBarry Smith   newC->localOffsets        = NULL;
1100298fd71SBarry Smith   newC->fillCount           = NULL;
1110298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1120298fd71SBarry Smith   newC->request             = NULL;
1130298fd71SBarry Smith   newC->status              = NULL;
114af33a6ddSJed Brown   PetscFunctionReturn(0);
115af33a6ddSJed Brown }
116af33a6ddSJed Brown 
117af33a6ddSJed Brown /*@C
118af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
119af33a6ddSJed Brown 
120af33a6ddSJed Brown    Logically Collective on Characteristic
121af33a6ddSJed Brown 
122af33a6ddSJed Brown    Input Parameters:
123af33a6ddSJed Brown +  c    - the method of characteristics context
124af33a6ddSJed Brown -  type - a known method
125af33a6ddSJed Brown 
126af33a6ddSJed Brown    Options Database Key:
127af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
128af33a6ddSJed Brown     of available methods
129af33a6ddSJed Brown 
130af33a6ddSJed Brown    Notes:
13130a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
132af33a6ddSJed Brown 
133af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
134af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
135af33a6ddSJed Brown   this routine.  Using the options database provides the user with
136af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
137af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
138af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
139af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
140af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
141af33a6ddSJed Brown   program, and the user's application is taking responsibility for
142af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
143af33a6ddSJed Brown   not for beginners.
144af33a6ddSJed Brown 
145af33a6ddSJed Brown   Level: intermediate
146af33a6ddSJed Brown 
147db781477SPatrick Sanan .seealso: `CharacteristicType`
148af33a6ddSJed Brown 
149af33a6ddSJed Brown @*/
15019fd82e9SBarry Smith PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
151af33a6ddSJed Brown {
152af33a6ddSJed Brown   PetscBool      match;
1535f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
154af33a6ddSJed Brown 
155af33a6ddSJed Brown   PetscFunctionBegin;
156af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
157af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
158af33a6ddSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) c, type, &match));
160af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
161af33a6ddSJed Brown 
162af33a6ddSJed Brown   if (c->data) {
163af33a6ddSJed Brown     /* destroy the old private Characteristic context */
164*dbbe0bcdSBarry Smith     PetscUseTypeMethod(c,destroy);
1650298fd71SBarry Smith     c->ops->destroy = NULL;
1662120983fSLisandro Dalcin     c->data         = NULL;
167af33a6ddSJed Brown   }
168af33a6ddSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(CharacteristicList,type,&r));
1703c633725SBarry Smith   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
171af33a6ddSJed Brown   c->setupcalled = 0;
1729566063dSJacob Faibussowitsch   PetscCall((*r)(c));
1739566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject) c, type));
174af33a6ddSJed Brown   PetscFunctionReturn(0);
175af33a6ddSJed Brown }
176af33a6ddSJed Brown 
177af33a6ddSJed Brown /*@
178af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
179af33a6ddSJed Brown    later use of an iterative solver.
180af33a6ddSJed Brown 
181af33a6ddSJed Brown    Collective on Characteristic
182af33a6ddSJed Brown 
183af33a6ddSJed Brown    Input Parameter:
184af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
185af33a6ddSJed Brown 
186af33a6ddSJed Brown    Level: developer
187af33a6ddSJed Brown 
188db781477SPatrick Sanan .seealso: `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
189af33a6ddSJed Brown @*/
190af33a6ddSJed Brown PetscErrorCode CharacteristicSetUp(Characteristic c)
191af33a6ddSJed Brown {
192af33a6ddSJed Brown   PetscFunctionBegin;
193af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
194af33a6ddSJed Brown 
195af33a6ddSJed Brown   if (!((PetscObject)c)->type_name) {
1969566063dSJacob Faibussowitsch     PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
197af33a6ddSJed Brown   }
198af33a6ddSJed Brown 
199af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
200af33a6ddSJed Brown 
2019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL));
202af33a6ddSJed Brown   if (!c->setupcalled) {
203*dbbe0bcdSBarry Smith     PetscUseTypeMethod(c,setup);
204af33a6ddSJed Brown   }
2059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp,c,NULL,NULL,NULL));
206af33a6ddSJed Brown   c->setupcalled = 2;
207af33a6ddSJed Brown   PetscFunctionReturn(0);
208af33a6ddSJed Brown }
209af33a6ddSJed Brown 
210af33a6ddSJed Brown /*@C
2111c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2121c84c290SBarry Smith 
2131c84c290SBarry Smith    Not Collective
2141c84c290SBarry Smith 
2151c84c290SBarry Smith    Input Parameters:
2161c84c290SBarry Smith +  name_solver - name of a new user-defined solver
2171c84c290SBarry Smith -  routine_create - routine to create method context
2181c84c290SBarry Smith 
2191c84c290SBarry Smith   Sample usage:
2201c84c290SBarry Smith .vb
221bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2221c84c290SBarry Smith .ve
2231c84c290SBarry Smith 
2241c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2251c84c290SBarry Smith .vb
2261c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2271c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2281c84c290SBarry Smith .ve
2291c84c290SBarry Smith    or at runtime via the option
2301c84c290SBarry Smith .vb
2311c84c290SBarry Smith     -characteristic_type my_char
2321c84c290SBarry Smith .ve
2331c84c290SBarry Smith 
2341c84c290SBarry Smith    Notes:
2351c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2361c84c290SBarry Smith 
237db781477SPatrick Sanan .seealso: `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
238af33a6ddSJed Brown 
239af33a6ddSJed Brown   Level: advanced
240af33a6ddSJed Brown @*/
241bdf89e91SBarry Smith PetscErrorCode CharacteristicRegister(const char sname[],PetscErrorCode (*function)(Characteristic))
242af33a6ddSJed Brown {
243af33a6ddSJed Brown   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
2459566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&CharacteristicList,sname,function));
246af33a6ddSJed Brown   PetscFunctionReturn(0);
247af33a6ddSJed Brown }
248af33a6ddSJed Brown 
249af33a6ddSJed 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)
250af33a6ddSJed Brown {
251af33a6ddSJed Brown   PetscFunctionBegin;
252af33a6ddSJed Brown   c->velocityDA      = da;
253af33a6ddSJed Brown   c->velocity        = v;
254af33a6ddSJed Brown   c->velocityOld     = vOld;
255af33a6ddSJed Brown   c->numVelocityComp = numComponents;
256af33a6ddSJed Brown   c->velocityComp    = components;
257af33a6ddSJed Brown   c->velocityInterp  = interp;
258af33a6ddSJed Brown   c->velocityCtx     = ctx;
259af33a6ddSJed Brown   PetscFunctionReturn(0);
260af33a6ddSJed Brown }
261af33a6ddSJed Brown 
262af33a6ddSJed 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)
263af33a6ddSJed Brown {
264af33a6ddSJed Brown   PetscFunctionBegin;
265af33a6ddSJed Brown   c->velocityDA          = da;
266af33a6ddSJed Brown   c->velocity            = v;
267af33a6ddSJed Brown   c->velocityOld         = vOld;
268af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
269af33a6ddSJed Brown   c->velocityComp        = components;
270af33a6ddSJed Brown   c->velocityInterpLocal = interp;
271af33a6ddSJed Brown   c->velocityCtx         = ctx;
272af33a6ddSJed Brown   PetscFunctionReturn(0);
273af33a6ddSJed Brown }
274af33a6ddSJed Brown 
275af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void*), void *ctx)
276af33a6ddSJed Brown {
277af33a6ddSJed Brown   PetscFunctionBegin;
278af33a6ddSJed Brown #if 0
2793c633725SBarry Smith   PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
280af33a6ddSJed Brown #endif
281af33a6ddSJed Brown   c->fieldDA      = da;
282af33a6ddSJed Brown   c->field        = v;
283af33a6ddSJed Brown   c->numFieldComp = numComponents;
284af33a6ddSJed Brown   c->fieldComp    = components;
285af33a6ddSJed Brown   c->fieldInterp  = interp;
286af33a6ddSJed Brown   c->fieldCtx     = ctx;
287af33a6ddSJed Brown   PetscFunctionReturn(0);
288af33a6ddSJed Brown }
289af33a6ddSJed Brown 
290af33a6ddSJed Brown PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void*, PetscReal[], PetscInt, PetscInt[], PetscScalar [], void*), void *ctx)
291af33a6ddSJed Brown {
292af33a6ddSJed Brown   PetscFunctionBegin;
293af33a6ddSJed Brown #if 0
2943c633725SBarry Smith   PetscCheck(numComponents <= 2,PETSC_COMM_SELF,PETSC_ERR_SUP, "Fields with more than 2 components are not supported. Send mail to petsc-maint@mcs.anl.gov.");
295af33a6ddSJed Brown #endif
296af33a6ddSJed Brown   c->fieldDA          = da;
297af33a6ddSJed Brown   c->field            = v;
298af33a6ddSJed Brown   c->numFieldComp     = numComponents;
299af33a6ddSJed Brown   c->fieldComp        = components;
300af33a6ddSJed Brown   c->fieldInterpLocal = interp;
301af33a6ddSJed Brown   c->fieldCtx         = ctx;
302af33a6ddSJed Brown   PetscFunctionReturn(0);
303af33a6ddSJed Brown }
304af33a6ddSJed Brown 
305af33a6ddSJed Brown PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
306af33a6ddSJed Brown {
307af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
308af33a6ddSJed Brown   DM                      da = c->velocityDA;
309af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
310af33a6ddSJed Brown   Vec                     fieldLocal;
311af33a6ddSJed Brown   DMDALocalInfo           info;
312af33a6ddSJed Brown   PetscScalar             **solArray;
313af33a6ddSJed Brown   void                    *velocityArray;
314af33a6ddSJed Brown   void                    *velocityArrayOld;
315af33a6ddSJed Brown   void                    *fieldArray;
3161cc8b206SBarry Smith   PetscScalar             *interpIndices;
3171cc8b206SBarry Smith   PetscScalar             *velocityValues, *velocityValuesOld;
3181cc8b206SBarry Smith   PetscScalar             *fieldValues;
319af33a6ddSJed Brown   PetscMPIInt             rank;
320af33a6ddSJed Brown   PetscInt                dim;
321af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
322af33a6ddSJed Brown   PetscInt                dof;
323af33a6ddSJed Brown   PetscInt                gx, gy;
324af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
325af33a6ddSJed Brown 
326af33a6ddSJed Brown   PetscFunctionBegin;
327af33a6ddSJed Brown   c->queueSize = 0;
3289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3299566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighborsRank(da, neighbors));
3309566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3319566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetUp(c));
332af33a6ddSJed Brown   /* global and local grid info */
3339566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3349566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
335af33a6ddSJed Brown   is   = info.xs;          ie   = info.xs+info.xm;
336af33a6ddSJed Brown   js   = info.ys;          je   = info.ys+info.ym;
337af33a6ddSJed Brown   /* Allocation */
3389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim,                &interpIndices));
3399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numFieldComp,    &fieldValues));
3429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL));
343af33a6ddSJed Brown 
344af33a6ddSJed Brown   /* -----------------------------------------------------------------------
345af33a6ddSJed Brown      PART 1, AT t-dt/2
346af33a6ddSJed Brown      -----------------------------------------------------------------------*/
3479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL));
348af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
349af33a6ddSJed Brown   if (c->velocityInterpLocal) {
3509566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3519566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3529566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3539566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3549566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3559566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3569566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal,    &velocityArray));
3579566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
358af33a6ddSJed Brown   }
3599566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
360af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
361af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
362af33a6ddSJed Brown       interpIndices[0] = Qi.i;
363af33a6ddSJed Brown       interpIndices[1] = Qi.j;
3649566063dSJacob Faibussowitsch       if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3659566063dSJacob Faibussowitsch       else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
366af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0]*dt/2.0;
367af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1]*dt/2.0;
368af33a6ddSJed Brown 
369af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
370af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
371af33a6ddSJed Brown 
372af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
3739566063dSJacob Faibussowitsch       PetscCall(DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y)));
3749566063dSJacob Faibussowitsch       PetscCall(CharacteristicAddPoint(c, &Qi));
375af33a6ddSJed Brown     }
376af33a6ddSJed Brown   }
3779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup,NULL,NULL,NULL,NULL));
378af33a6ddSJed Brown 
3799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
3809566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
3819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
382af33a6ddSJed Brown 
3839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL));
384af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3859566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
386af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
387af33a6ddSJed Brown     Qi = c->queue[n];
388af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
389af33a6ddSJed Brown       interpIndices[0] = Qi.x;
390af33a6ddSJed Brown       interpIndices[1] = Qi.y;
391af33a6ddSJed Brown       if (c->velocityInterpLocal) {
3929566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3939566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
394af33a6ddSJed Brown       } else {
3959566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3969566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
397af33a6ddSJed Brown       }
398af33a6ddSJed Brown       Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
399af33a6ddSJed Brown       Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
400af33a6ddSJed Brown     }
401af33a6ddSJed Brown     c->queue[n] = Qi;
402af33a6ddSJed Brown   }
4039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal,NULL,NULL,NULL,NULL));
404af33a6ddSJed Brown 
4059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
4069566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
408af33a6ddSJed Brown 
409af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
4109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL));
41163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
412af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
413af33a6ddSJed Brown     Qi = c->queueRemote[n];
414af33a6ddSJed Brown     interpIndices[0] = Qi.x;
415af33a6ddSJed Brown     interpIndices[1] = Qi.y;
416af33a6ddSJed Brown     if (c->velocityInterpLocal) {
4179566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArray,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx));
4189566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
419af33a6ddSJed Brown     } else {
4209566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocity,    interpIndices, c->numVelocityComp, c->velocityComp, velocityValues,    c->velocityCtx));
4219566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
422af33a6ddSJed Brown     }
423af33a6ddSJed Brown     Qi.x = 0.5*(velocityValues[0] + velocityValuesOld[0]);
424af33a6ddSJed Brown     Qi.y = 0.5*(velocityValues[1] + velocityValuesOld[1]);
425af33a6ddSJed Brown     c->queueRemote[n] = Qi;
426af33a6ddSJed Brown   }
4279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote,NULL,NULL,NULL,NULL));
4289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
4299566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4309566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
431af33a6ddSJed Brown   if (c->velocityInterpLocal) {
4329566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal,    &velocityArray));
4339566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4349566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4359566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
436af33a6ddSJed Brown   }
4379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange,NULL,NULL,NULL,NULL));
438af33a6ddSJed Brown 
439af33a6ddSJed Brown   /* -----------------------------------------------------------------------
440af33a6ddSJed Brown      PART 2, AT t-dt
441af33a6ddSJed Brown      -----------------------------------------------------------------------*/
442af33a6ddSJed Brown 
443af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
4449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
4459566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
446af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
447af33a6ddSJed Brown     Qi   = c->queue[n];
448af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x*dt;
449af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y*dt;
450af33a6ddSJed Brown 
451af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
452af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
453af33a6ddSJed Brown 
454af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
4559566063dSJacob Faibussowitsch     PetscCall(DMDAMapCoordsToPeriodicDomain(da,&(Qi.x),&(Qi.y)));
456af33a6ddSJed Brown 
457af33a6ddSJed Brown     c->queue[n] = Qi;
458af33a6ddSJed Brown   }
4599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
460af33a6ddSJed Brown 
4619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
4629566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
4639566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
464af33a6ddSJed Brown 
465af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
467af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4689566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4699566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4709566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4719566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
472af33a6ddSJed Brown   }
4739566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
474af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
475af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
476af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
477af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
4789566063dSJacob Faibussowitsch       if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4799566063dSJacob Faibussowitsch       else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
480bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
481af33a6ddSJed Brown     }
482af33a6ddSJed Brown   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal,NULL,NULL,NULL,NULL));
484af33a6ddSJed Brown 
4859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
4869566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
488af33a6ddSJed Brown 
489af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
4909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL));
49163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
492af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
493af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
494af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
495af33a6ddSJed Brown 
496af33a6ddSJed Brown     /* for debugging purposes */
497af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
498af33a6ddSJed Brown       PetscScalar im = interpIndices[0]; PetscScalar jm = interpIndices[1];
499af33a6ddSJed Brown 
50063a3b9bcSJacob Faibussowitsch       PetscCheck((im >= (PetscScalar) is - 1.) && (im <= (PetscScalar) ie) && (jm >= (PetscScalar)  js - 1.) && (jm <= (PetscScalar) je),PETSC_COMM_SELF,PETSC_ERR_LIB, "Nonlocal point: (%g,%g)", (double)PetscAbsScalar(im), (double)PetscAbsScalar(jm));
501af33a6ddSJed Brown     }
502af33a6ddSJed Brown 
5039566063dSJacob Faibussowitsch     if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
5049566063dSJacob Faibussowitsch     else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
505bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
506af33a6ddSJed Brown   }
5079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote,NULL,NULL,NULL,NULL));
508af33a6ddSJed Brown 
5099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
5109566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
5119566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
512af33a6ddSJed Brown   if (c->fieldInterpLocal) {
5139566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
5149566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
515af33a6ddSJed Brown   }
5169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange,NULL,NULL,NULL,NULL));
517af33a6ddSJed Brown 
518af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
5199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL));
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(c->fieldDA,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
5219566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
522af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
523af33a6ddSJed Brown     Qi = c->queue[n];
524bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i*dof+c->fieldComp[comp]] = Qi.field[comp];
525af33a6ddSJed Brown   }
5269566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate,NULL,NULL,NULL,NULL));
5289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve,NULL,NULL,NULL,NULL));
529af33a6ddSJed Brown 
530af33a6ddSJed Brown   /* Cleanup */
5319566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpIndices));
5329566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValues));
5339566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValuesOld));
5349566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldValues));
535af33a6ddSJed Brown   PetscFunctionReturn(0);
536af33a6ddSJed Brown }
537af33a6ddSJed Brown 
538af33a6ddSJed Brown PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
539af33a6ddSJed Brown {
540af33a6ddSJed Brown   PetscFunctionBegin;
541af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5429566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->neighbors));
5439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5449566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
545af33a6ddSJed Brown   PetscFunctionReturn(0);
546af33a6ddSJed Brown }
547af33a6ddSJed Brown 
548af33a6ddSJed Brown PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
549af33a6ddSJed Brown {
550af33a6ddSJed Brown   PetscFunctionBegin;
55163a3b9bcSJacob Faibussowitsch   PetscCheck(c->queueSize < c->queueMax,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
552af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
553af33a6ddSJed Brown   PetscFunctionReturn(0);
554af33a6ddSJed Brown }
555af33a6ddSJed Brown 
556af33a6ddSJed Brown int CharacteristicSendCoordinatesBegin(Characteristic c)
557af33a6ddSJed Brown {
558af33a6ddSJed Brown   PetscMPIInt    rank, tag = 121;
559af33a6ddSJed Brown   PetscInt       i, n;
560af33a6ddSJed Brown 
561af33a6ddSJed Brown   PetscFunctionBegin;
5629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5639566063dSJacob Faibussowitsch   PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5649566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
565bbd56ea5SKarl Rupp   for (i = 0;  i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
566af33a6ddSJed Brown   c->fillCount[0] = 0;
567af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
5689566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
569af33a6ddSJed Brown   }
570af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
5719566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
572af33a6ddSJed Brown   }
5739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
574af33a6ddSJed Brown   /* Initialize the remote queue */
575af33a6ddSJed Brown   c->queueLocalMax  = c->localOffsets[0]  = 0;
576af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
577af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
578af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
579af33a6ddSJed Brown     c->queueRemoteMax  += c->fillCount[n];
580af33a6ddSJed Brown     c->localOffsets[n]  = c->queueLocalMax;
581af33a6ddSJed Brown     c->queueLocalMax   += c->needCount[n];
582af33a6ddSJed Brown   }
583af33a6ddSJed Brown   /* HACK BEGIN */
584bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
585af33a6ddSJed Brown   c->needCount[0] = 0;
586af33a6ddSJed Brown   /* HACK END */
587af33a6ddSJed Brown   if (c->queueRemoteMax) {
5889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5890298fd71SBarry Smith   } else c->queueRemote = NULL;
590af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
591af33a6ddSJed Brown 
592af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
593af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
59463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
5959566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
596af33a6ddSJed Brown   }
597af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
59863a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
5999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
600af33a6ddSJed Brown   }
601af33a6ddSJed Brown   PetscFunctionReturn(0);
602af33a6ddSJed Brown }
603af33a6ddSJed Brown 
604af33a6ddSJed Brown PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
605af33a6ddSJed Brown {
606af33a6ddSJed Brown #if 0
607af33a6ddSJed Brown   PetscMPIInt rank;
608af33a6ddSJed Brown   PetscInt    n;
609af33a6ddSJed Brown #endif
610af33a6ddSJed Brown 
611af33a6ddSJed Brown   PetscFunctionBegin;
6129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
613af33a6ddSJed Brown #if 0
6149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
615af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
6163c633725SBarry Smith     PetscCheck(c->neighbors[c->queueRemote[n].proc] != rank,PETSC_COMM_SELF,PETSC_ERR_PLIB, "This is messed up, n = %d proc = %d", n, c->queueRemote[n].proc);
617af33a6ddSJed Brown   }
618af33a6ddSJed Brown #endif
619af33a6ddSJed Brown   PetscFunctionReturn(0);
620af33a6ddSJed Brown }
621af33a6ddSJed Brown 
622af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
623af33a6ddSJed Brown {
624af33a6ddSJed Brown   PetscMPIInt    tag = 121;
625af33a6ddSJed Brown   PetscInt       n;
626af33a6ddSJed Brown 
627af33a6ddSJed Brown   PetscFunctionBegin;
628af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
629af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6309566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n-1])));
631af33a6ddSJed Brown   }
632af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
6339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
634af33a6ddSJed Brown   }
635af33a6ddSJed Brown   PetscFunctionReturn(0);
636af33a6ddSJed Brown }
637af33a6ddSJed Brown 
638af33a6ddSJed Brown PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
639af33a6ddSJed Brown {
640af33a6ddSJed Brown   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors-1, c->request, c->status));
642af33a6ddSJed Brown   /* Free queue of requests from other procs */
6439566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->queueRemote));
644af33a6ddSJed Brown   PetscFunctionReturn(0);
645af33a6ddSJed Brown }
646af33a6ddSJed Brown 
647af33a6ddSJed Brown /*---------------------------------------------------------------------*/
648af33a6ddSJed Brown /*
649af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
650af33a6ddSJed Brown */
6510b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
652af33a6ddSJed Brown /*---------------------------------------------------------------------*/
653af33a6ddSJed Brown {
654af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6550b8f38c8SBarry Smith   PetscInt                n;
656af33a6ddSJed Brown 
6570b8f38c8SBarry Smith   PetscFunctionBegin;
658af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6599566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
660af33a6ddSJed Brown     for (n=0; n<size; n++) {
66163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"%" PetscInt_FMT " %d\n",n,queue[n].proc));
662af33a6ddSJed Brown     }
663af33a6ddSJed Brown   }
664af33a6ddSJed Brown 
665af33a6ddSJed Brown   /* SORTING PHASE */
6660b8f38c8SBarry Smith   for (n = (size / 2)-1; n >= 0; n--) {
6679566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, n, size-1)); /* Rich had size-1 here, Matt had size*/
6680b8f38c8SBarry Smith   }
669af33a6ddSJed Brown   for (n = size-1; n >= 1; n--) {
670af33a6ddSJed Brown     temp     = queue[0];
671af33a6ddSJed Brown     queue[0] = queue[n];
672af33a6ddSJed Brown     queue[n] = temp;
6739566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, 0, n-1));
674af33a6ddSJed Brown   }
675af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
6769566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Avter  Heap sort\n"));
677af33a6ddSJed Brown     for (n=0; n<size; n++) {
67863a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(NULL,"%" PetscInt_FMT " %d\n",n,queue[n].proc));
679af33a6ddSJed Brown     }
680af33a6ddSJed Brown   }
6810b8f38c8SBarry Smith   PetscFunctionReturn(0);
682af33a6ddSJed Brown }
683af33a6ddSJed Brown 
684af33a6ddSJed Brown /*---------------------------------------------------------------------*/
685af33a6ddSJed Brown /*
686af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
687af33a6ddSJed Brown */
6880b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
689af33a6ddSJed Brown /*---------------------------------------------------------------------*/
690af33a6ddSJed Brown {
691af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
692af33a6ddSJed Brown   PetscInt                maxChild;
693af33a6ddSJed Brown   CharacteristicPointDA2D temp;
694af33a6ddSJed Brown 
695af33a6ddSJed Brown   PetscFunctionBegin;
696af33a6ddSJed Brown   while ((root*2 <= bottom) && (!done)) {
697af33a6ddSJed Brown     if (root*2 == bottom) maxChild = root * 2;
698af33a6ddSJed Brown     else if (queue[root*2].proc > queue[root*2+1].proc) maxChild = root * 2;
699af33a6ddSJed Brown     else maxChild = root * 2 + 1;
700af33a6ddSJed Brown 
701af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
702af33a6ddSJed Brown       temp            = queue[root];
703af33a6ddSJed Brown       queue[root]     = queue[maxChild];
704af33a6ddSJed Brown       queue[maxChild] = temp;
705af33a6ddSJed Brown       root            = maxChild;
706db4deed7SKarl Rupp     } else done = PETSC_TRUE;
707af33a6ddSJed Brown   }
708af33a6ddSJed Brown   PetscFunctionReturn(0);
709af33a6ddSJed Brown }
710af33a6ddSJed Brown 
711af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
712af33a6ddSJed Brown PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
713af33a6ddSJed Brown {
714bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx, by;
715af33a6ddSJed Brown   PetscBool        IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
716af33a6ddSJed Brown   MPI_Comm         comm;
717af33a6ddSJed Brown   PetscMPIInt      rank;
718af33a6ddSJed Brown   PetscInt         **procs,pi,pj,pim,pip,pjm,pjp,PI,PJ;
719af33a6ddSJed Brown 
720af33a6ddSJed Brown   PetscFunctionBegin;
7219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) da, &comm));
7229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7239566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI,&PJ, NULL, NULL, NULL, &bx, &by,NULL, NULL));
724af33a6ddSJed Brown 
725bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
726bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
727af33a6ddSJed Brown 
728af33a6ddSJed Brown   neighbors[0] = rank;
729af33a6ddSJed Brown   rank = 0;
7309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PJ,&procs));
731af33a6ddSJed Brown   for (pj=0; pj<PJ; pj++) {
7329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PI,&(procs[pj])));
733af33a6ddSJed Brown     for (pi=0; pi<PI; pi++) {
734af33a6ddSJed Brown       procs[pj][pi] = rank;
735af33a6ddSJed Brown       rank++;
736af33a6ddSJed Brown     }
737af33a6ddSJed Brown   }
738af33a6ddSJed Brown 
739af33a6ddSJed Brown   pi  = neighbors[0] % PI;
740af33a6ddSJed Brown   pj  = neighbors[0] / PI;
741af33a6ddSJed Brown   pim = pi-1;  if (pim<0) pim=PI-1;
742af33a6ddSJed Brown   pip = (pi+1)%PI;
743af33a6ddSJed Brown   pjm = pj-1;  if (pjm<0) pjm=PJ-1;
744af33a6ddSJed Brown   pjp = (pj+1)%PJ;
745af33a6ddSJed Brown 
746af33a6ddSJed Brown   neighbors[1] = procs[pj] [pim];
747af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
748af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
749af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
750af33a6ddSJed Brown   neighbors[5] = procs[pj] [pip];
751af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
752af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
753af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
754af33a6ddSJed Brown 
755af33a6ddSJed Brown   if (!IPeriodic) {
756af33a6ddSJed Brown     if (pi==0)    neighbors[1]=neighbors[2]=neighbors[8]=neighbors[0];
757af33a6ddSJed Brown     if (pi==PI-1) neighbors[4]=neighbors[5]=neighbors[6]=neighbors[0];
758af33a6ddSJed Brown   }
759af33a6ddSJed Brown 
760af33a6ddSJed Brown   if (!JPeriodic) {
761af33a6ddSJed Brown     if (pj==0)    neighbors[6]=neighbors[7]=neighbors[8]=neighbors[0];
762af33a6ddSJed Brown     if (pj==PJ-1) neighbors[2]=neighbors[3]=neighbors[4]=neighbors[0];
763af33a6ddSJed Brown   }
764af33a6ddSJed Brown 
765af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
7669566063dSJacob Faibussowitsch     PetscCall(PetscFree(procs[pj]));
767af33a6ddSJed Brown   }
7689566063dSJacob Faibussowitsch   PetscCall(PetscFree(procs));
769af33a6ddSJed Brown   PetscFunctionReturn(0);
770af33a6ddSJed Brown }
771af33a6ddSJed Brown 
772af33a6ddSJed Brown /*
773af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
774af33a6ddSJed Brown     2 | 3 | 4
775af33a6ddSJed Brown     __|___|__
776af33a6ddSJed Brown     1 | 0 | 5
777af33a6ddSJed Brown     __|___|__
778af33a6ddSJed Brown     8 | 7 | 6
779af33a6ddSJed Brown       |   |
780af33a6ddSJed Brown */
7811cc8b206SBarry Smith PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
782af33a6ddSJed Brown {
783af33a6ddSJed Brown   DMDALocalInfo  info;
7841cc8b206SBarry Smith   PetscReal      is,ie,js,je;
785af33a6ddSJed Brown 
7869566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
7871cc8b206SBarry Smith   is   = (PetscReal) info.xs - 0.5; ie = (PetscReal) info.xs + info.xm - 0.5;
7881cc8b206SBarry Smith   js   = (PetscReal) info.ys - 0.5; je = (PetscReal) info.ys + info.ym - 0.5;
789af33a6ddSJed Brown 
790af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
791bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
792bbd56ea5SKarl Rupp     else if (jr < js)         return 7;
793bbd56ea5SKarl Rupp     else                      return 3;
794af33a6ddSJed Brown   } else if (ir < is) {     /* left column */
795bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
796bbd56ea5SKarl Rupp     else if (jr < js)         return 8;
797bbd56ea5SKarl Rupp     else                      return 2;
798af33a6ddSJed Brown   } else {                  /* right column */
799bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
800bbd56ea5SKarl Rupp     else if (jr < js)         return 6;
801bbd56ea5SKarl Rupp     else                      return 4;
802af33a6ddSJed Brown   }
803af33a6ddSJed Brown }
804