xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
23*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
24*d71ae5a4SJacob Faibussowitsch {
25af33a6ddSJed Brown   PetscBool iascii;
26af33a6ddSJed Brown 
27af33a6ddSJed Brown   PetscFunctionBegin;
28af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
2948a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
30af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
31af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
32af33a6ddSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
34ad540459SPierre Jolivet   if (!iascii) PetscTryTypeMethod(c, view, viewer);
35af33a6ddSJed Brown   PetscFunctionReturn(0);
36af33a6ddSJed Brown }
37af33a6ddSJed Brown 
38*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicDestroy(Characteristic *c)
39*d71ae5a4SJacob Faibussowitsch {
40af33a6ddSJed Brown   PetscFunctionBegin;
416bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
426bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
436bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(0);
44af33a6ddSJed Brown 
4548a46eb9SPierre Jolivet   if ((*c)->ops->destroy) PetscCall((*(*c)->ops->destroy)((*c)));
469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&(*c)->itemType));
479566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queue));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueLocal));
499566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueRemote));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->neighbors));
519566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->needCount));
529566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->localOffsets));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->fillCount));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->remoteOffsets));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->request));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->status));
579566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
58af33a6ddSJed Brown   PetscFunctionReturn(0);
59af33a6ddSJed Brown }
60af33a6ddSJed Brown 
61*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
62*d71ae5a4SJacob Faibussowitsch {
63af33a6ddSJed Brown   Characteristic newC;
64af33a6ddSJed Brown 
65af33a6ddSJed Brown   PetscFunctionBegin;
66af33a6ddSJed Brown   PetscValidPointer(c, 2);
670298fd71SBarry Smith   *c = NULL;
689566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
69af33a6ddSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
71af33a6ddSJed Brown   *c = newC;
72af33a6ddSJed Brown 
73af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
74af33a6ddSJed Brown   newC->numIds              = 0;
750298fd71SBarry Smith   newC->velocityDA          = NULL;
760298fd71SBarry Smith   newC->velocity            = NULL;
770298fd71SBarry Smith   newC->velocityOld         = NULL;
78af33a6ddSJed Brown   newC->numVelocityComp     = 0;
790298fd71SBarry Smith   newC->velocityComp        = NULL;
800298fd71SBarry Smith   newC->velocityInterp      = NULL;
810298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
820298fd71SBarry Smith   newC->velocityCtx         = NULL;
830298fd71SBarry Smith   newC->fieldDA             = NULL;
840298fd71SBarry Smith   newC->field               = NULL;
85af33a6ddSJed Brown   newC->numFieldComp        = 0;
860298fd71SBarry Smith   newC->fieldComp           = NULL;
870298fd71SBarry Smith   newC->fieldInterp         = NULL;
880298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
890298fd71SBarry Smith   newC->fieldCtx            = NULL;
900298fd71SBarry Smith   newC->itemType            = 0;
910298fd71SBarry Smith   newC->queue               = NULL;
92af33a6ddSJed Brown   newC->queueSize           = 0;
93af33a6ddSJed Brown   newC->queueMax            = 0;
940298fd71SBarry Smith   newC->queueLocal          = NULL;
95af33a6ddSJed Brown   newC->queueLocalSize      = 0;
96af33a6ddSJed Brown   newC->queueLocalMax       = 0;
970298fd71SBarry Smith   newC->queueRemote         = NULL;
98af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
99af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
100af33a6ddSJed Brown   newC->numNeighbors        = 0;
1010298fd71SBarry Smith   newC->neighbors           = NULL;
1020298fd71SBarry Smith   newC->needCount           = NULL;
1030298fd71SBarry Smith   newC->localOffsets        = NULL;
1040298fd71SBarry Smith   newC->fillCount           = NULL;
1050298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1060298fd71SBarry Smith   newC->request             = NULL;
1070298fd71SBarry Smith   newC->status              = NULL;
108af33a6ddSJed Brown   PetscFunctionReturn(0);
109af33a6ddSJed Brown }
110af33a6ddSJed Brown 
111af33a6ddSJed Brown /*@C
112af33a6ddSJed Brown    CharacteristicSetType - Builds Characteristic for a particular solver.
113af33a6ddSJed Brown 
114af33a6ddSJed Brown    Logically Collective on Characteristic
115af33a6ddSJed Brown 
116af33a6ddSJed Brown    Input Parameters:
117af33a6ddSJed Brown +  c    - the method of characteristics context
118af33a6ddSJed Brown -  type - a known method
119af33a6ddSJed Brown 
120af33a6ddSJed Brown    Options Database Key:
121af33a6ddSJed Brown .  -characteristic_type <method> - Sets the method; use -help for a list
122af33a6ddSJed Brown     of available methods
123af33a6ddSJed Brown 
124af33a6ddSJed Brown    Notes:
12530a76a96SBarry Smith    See "include/petsccharacteristic.h" for available methods
126af33a6ddSJed Brown 
127af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
128af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
129af33a6ddSJed Brown   this routine.  Using the options database provides the user with
130af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
131af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
132af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
133af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
134af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
135af33a6ddSJed Brown   program, and the user's application is taking responsibility for
136af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
137af33a6ddSJed Brown   not for beginners.
138af33a6ddSJed Brown 
139af33a6ddSJed Brown   Level: intermediate
140af33a6ddSJed Brown 
141db781477SPatrick Sanan .seealso: `CharacteristicType`
142af33a6ddSJed Brown 
143af33a6ddSJed Brown @*/
144*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
145*d71ae5a4SJacob Faibussowitsch {
146af33a6ddSJed Brown   PetscBool match;
1475f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
148af33a6ddSJed Brown 
149af33a6ddSJed Brown   PetscFunctionBegin;
150af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
151af33a6ddSJed Brown   PetscValidCharPointer(type, 2);
152af33a6ddSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match));
154af33a6ddSJed Brown   if (match) PetscFunctionReturn(0);
155af33a6ddSJed Brown 
156af33a6ddSJed Brown   if (c->data) {
157af33a6ddSJed Brown     /* destroy the old private Characteristic context */
158dbbe0bcdSBarry Smith     PetscUseTypeMethod(c, destroy);
1590298fd71SBarry Smith     c->ops->destroy = NULL;
1602120983fSLisandro Dalcin     c->data         = NULL;
161af33a6ddSJed Brown   }
162af33a6ddSJed Brown 
1639566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(CharacteristicList, type, &r));
1643c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
165af33a6ddSJed Brown   c->setupcalled = 0;
1669566063dSJacob Faibussowitsch   PetscCall((*r)(c));
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)c, type));
168af33a6ddSJed Brown   PetscFunctionReturn(0);
169af33a6ddSJed Brown }
170af33a6ddSJed Brown 
171af33a6ddSJed Brown /*@
172af33a6ddSJed Brown    CharacteristicSetUp - Sets up the internal data structures for the
173af33a6ddSJed Brown    later use of an iterative solver.
174af33a6ddSJed Brown 
175af33a6ddSJed Brown    Collective on Characteristic
176af33a6ddSJed Brown 
177af33a6ddSJed Brown    Input Parameter:
178af33a6ddSJed Brown .  ksp   - iterative context obtained from CharacteristicCreate()
179af33a6ddSJed Brown 
180af33a6ddSJed Brown    Level: developer
181af33a6ddSJed Brown 
182db781477SPatrick Sanan .seealso: `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
183af33a6ddSJed Brown @*/
184*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetUp(Characteristic c)
185*d71ae5a4SJacob Faibussowitsch {
186af33a6ddSJed Brown   PetscFunctionBegin;
187af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
188af33a6ddSJed Brown 
18948a46eb9SPierre Jolivet   if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
190af33a6ddSJed Brown 
191af33a6ddSJed Brown   if (c->setupcalled == 2) PetscFunctionReturn(0);
192af33a6ddSJed Brown 
1939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
194ad540459SPierre Jolivet   if (!c->setupcalled) PetscUseTypeMethod(c, setup);
1959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
196af33a6ddSJed Brown   c->setupcalled = 2;
197af33a6ddSJed Brown   PetscFunctionReturn(0);
198af33a6ddSJed Brown }
199af33a6ddSJed Brown 
200af33a6ddSJed Brown /*@C
2011c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2021c84c290SBarry Smith 
2031c84c290SBarry Smith    Not Collective
2041c84c290SBarry Smith 
2051c84c290SBarry Smith    Input Parameters:
2061c84c290SBarry Smith +  name_solver - name of a new user-defined solver
2071c84c290SBarry Smith -  routine_create - routine to create method context
2081c84c290SBarry Smith 
2091c84c290SBarry Smith   Sample usage:
2101c84c290SBarry Smith .vb
211bdf89e91SBarry Smith     CharacteristicRegister("my_char", MyCharCreate);
2121c84c290SBarry Smith .ve
2131c84c290SBarry Smith 
2141c84c290SBarry Smith   Then, your Characteristic type can be chosen with the procedural interface via
2151c84c290SBarry Smith .vb
2161c84c290SBarry Smith     CharacteristicCreate(MPI_Comm, Characteristic* &char);
2171c84c290SBarry Smith     CharacteristicSetType(char,"my_char");
2181c84c290SBarry Smith .ve
2191c84c290SBarry Smith    or at runtime via the option
2201c84c290SBarry Smith .vb
2211c84c290SBarry Smith     -characteristic_type my_char
2221c84c290SBarry Smith .ve
2231c84c290SBarry Smith 
2241c84c290SBarry Smith    Notes:
2251c84c290SBarry Smith    CharacteristicRegister() may be called multiple times to add several user-defined solvers.
2261c84c290SBarry Smith 
227db781477SPatrick Sanan .seealso: `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
228af33a6ddSJed Brown 
229af33a6ddSJed Brown   Level: advanced
230af33a6ddSJed Brown @*/
231*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic))
232*d71ae5a4SJacob Faibussowitsch {
233af33a6ddSJed Brown   PetscFunctionBegin;
2349566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
2359566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function));
236af33a6ddSJed Brown   PetscFunctionReturn(0);
237af33a6ddSJed Brown }
238af33a6ddSJed Brown 
239*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetVelocityInterpolation(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
240*d71ae5a4SJacob Faibussowitsch {
241af33a6ddSJed Brown   PetscFunctionBegin;
242af33a6ddSJed Brown   c->velocityDA      = da;
243af33a6ddSJed Brown   c->velocity        = v;
244af33a6ddSJed Brown   c->velocityOld     = vOld;
245af33a6ddSJed Brown   c->numVelocityComp = numComponents;
246af33a6ddSJed Brown   c->velocityComp    = components;
247af33a6ddSJed Brown   c->velocityInterp  = interp;
248af33a6ddSJed Brown   c->velocityCtx     = ctx;
249af33a6ddSJed Brown   PetscFunctionReturn(0);
250af33a6ddSJed Brown }
251af33a6ddSJed Brown 
252*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetVelocityInterpolationLocal(Characteristic c, DM da, Vec v, Vec vOld, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
253*d71ae5a4SJacob Faibussowitsch {
254af33a6ddSJed Brown   PetscFunctionBegin;
255af33a6ddSJed Brown   c->velocityDA          = da;
256af33a6ddSJed Brown   c->velocity            = v;
257af33a6ddSJed Brown   c->velocityOld         = vOld;
258af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
259af33a6ddSJed Brown   c->velocityComp        = components;
260af33a6ddSJed Brown   c->velocityInterpLocal = interp;
261af33a6ddSJed Brown   c->velocityCtx         = ctx;
262af33a6ddSJed Brown   PetscFunctionReturn(0);
263af33a6ddSJed Brown }
264af33a6ddSJed Brown 
265*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
266*d71ae5a4SJacob Faibussowitsch {
267af33a6ddSJed Brown   PetscFunctionBegin;
268af33a6ddSJed Brown #if 0
2693c633725SBarry 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.");
270af33a6ddSJed Brown #endif
271af33a6ddSJed Brown   c->fieldDA      = da;
272af33a6ddSJed Brown   c->field        = v;
273af33a6ddSJed Brown   c->numFieldComp = numComponents;
274af33a6ddSJed Brown   c->fieldComp    = components;
275af33a6ddSJed Brown   c->fieldInterp  = interp;
276af33a6ddSJed Brown   c->fieldCtx     = ctx;
277af33a6ddSJed Brown   PetscFunctionReturn(0);
278af33a6ddSJed Brown }
279af33a6ddSJed Brown 
280*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
281*d71ae5a4SJacob Faibussowitsch {
282af33a6ddSJed Brown   PetscFunctionBegin;
283af33a6ddSJed Brown #if 0
2843c633725SBarry 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.");
285af33a6ddSJed Brown #endif
286af33a6ddSJed Brown   c->fieldDA          = da;
287af33a6ddSJed Brown   c->field            = v;
288af33a6ddSJed Brown   c->numFieldComp     = numComponents;
289af33a6ddSJed Brown   c->fieldComp        = components;
290af33a6ddSJed Brown   c->fieldInterpLocal = interp;
291af33a6ddSJed Brown   c->fieldCtx         = ctx;
292af33a6ddSJed Brown   PetscFunctionReturn(0);
293af33a6ddSJed Brown }
294af33a6ddSJed Brown 
295*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
296*d71ae5a4SJacob Faibussowitsch {
297af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
298af33a6ddSJed Brown   DM                      da = c->velocityDA;
299af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
300af33a6ddSJed Brown   Vec                     fieldLocal;
301af33a6ddSJed Brown   DMDALocalInfo           info;
302af33a6ddSJed Brown   PetscScalar           **solArray;
303af33a6ddSJed Brown   void                   *velocityArray;
304af33a6ddSJed Brown   void                   *velocityArrayOld;
305af33a6ddSJed Brown   void                   *fieldArray;
3061cc8b206SBarry Smith   PetscScalar            *interpIndices;
3071cc8b206SBarry Smith   PetscScalar            *velocityValues, *velocityValuesOld;
3081cc8b206SBarry Smith   PetscScalar            *fieldValues;
309af33a6ddSJed Brown   PetscMPIInt             rank;
310af33a6ddSJed Brown   PetscInt                dim;
311af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
312af33a6ddSJed Brown   PetscInt                dof;
313af33a6ddSJed Brown   PetscInt                gx, gy;
314af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
315af33a6ddSJed Brown 
316af33a6ddSJed Brown   PetscFunctionBegin;
317af33a6ddSJed Brown   c->queueSize = 0;
3189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3199566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighborsRank(da, neighbors));
3209566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3219566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetUp(c));
322af33a6ddSJed Brown   /* global and local grid info */
3239566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3249566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
3259371c9d4SSatish Balay   is = info.xs;
3269371c9d4SSatish Balay   ie = info.xs + info.xm;
3279371c9d4SSatish Balay   js = info.ys;
3289371c9d4SSatish Balay   je = info.ys + info.ym;
329af33a6ddSJed Brown   /* Allocation */
3309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &interpIndices));
3319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues));
3349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
335af33a6ddSJed Brown 
336af33a6ddSJed Brown   /* -----------------------------------------------------------------------
337af33a6ddSJed Brown      PART 1, AT t-dt/2
338af33a6ddSJed Brown      -----------------------------------------------------------------------*/
3399566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
340af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
341af33a6ddSJed Brown   if (c->velocityInterpLocal) {
3429566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3439566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3449566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3459566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3469566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3479566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3489566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray));
3499566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
350af33a6ddSJed Brown   }
3519566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
352af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
353af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
354af33a6ddSJed Brown       interpIndices[0] = Qi.i;
355af33a6ddSJed Brown       interpIndices[1] = Qi.j;
3569566063dSJacob Faibussowitsch       if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3579566063dSJacob Faibussowitsch       else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
358af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0] * dt / 2.0;
359af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1] * dt / 2.0;
360af33a6ddSJed Brown 
361af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
362af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
363af33a6ddSJed Brown 
364af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
3659566063dSJacob Faibussowitsch       PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
3669566063dSJacob Faibussowitsch       PetscCall(CharacteristicAddPoint(c, &Qi));
367af33a6ddSJed Brown     }
368af33a6ddSJed Brown   }
3699566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
370af33a6ddSJed Brown 
3719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3729566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
3739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
374af33a6ddSJed Brown 
3759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
376af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3779566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
378af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
379af33a6ddSJed Brown     Qi = c->queue[n];
380af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
381af33a6ddSJed Brown       interpIndices[0] = Qi.x;
382af33a6ddSJed Brown       interpIndices[1] = Qi.y;
383af33a6ddSJed Brown       if (c->velocityInterpLocal) {
3849566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3859566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
386af33a6ddSJed Brown       } else {
3879566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3889566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
389af33a6ddSJed Brown       }
390af33a6ddSJed Brown       Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
391af33a6ddSJed Brown       Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
392af33a6ddSJed Brown     }
393af33a6ddSJed Brown     c->queue[n] = Qi;
394af33a6ddSJed Brown   }
3959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
396af33a6ddSJed Brown 
3979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3989566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
3999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
400af33a6ddSJed Brown 
401af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
4029566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
40363a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
404af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
405af33a6ddSJed Brown     Qi               = c->queueRemote[n];
406af33a6ddSJed Brown     interpIndices[0] = Qi.x;
407af33a6ddSJed Brown     interpIndices[1] = Qi.y;
408af33a6ddSJed Brown     if (c->velocityInterpLocal) {
4099566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4109566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
411af33a6ddSJed Brown     } else {
4129566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4139566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
414af33a6ddSJed Brown     }
415af33a6ddSJed Brown     Qi.x              = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
416af33a6ddSJed Brown     Qi.y              = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
417af33a6ddSJed Brown     c->queueRemote[n] = Qi;
418af33a6ddSJed Brown   }
4199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
4209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4219566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4229566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
423af33a6ddSJed Brown   if (c->velocityInterpLocal) {
4249566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray));
4259566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4269566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4279566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
428af33a6ddSJed Brown   }
4299566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
430af33a6ddSJed Brown 
431af33a6ddSJed Brown   /* -----------------------------------------------------------------------
432af33a6ddSJed Brown      PART 2, AT t-dt
433af33a6ddSJed Brown      -----------------------------------------------------------------------*/
434af33a6ddSJed Brown 
435af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
4369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
4379566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
438af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
439af33a6ddSJed Brown     Qi   = c->queue[n];
440af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x * dt;
441af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y * dt;
442af33a6ddSJed Brown 
443af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
444af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
445af33a6ddSJed Brown 
446af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
4479566063dSJacob Faibussowitsch     PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
448af33a6ddSJed Brown 
449af33a6ddSJed Brown     c->queue[n] = Qi;
450af33a6ddSJed Brown   }
4519566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
452af33a6ddSJed Brown 
4539566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4549566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
4559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
456af33a6ddSJed Brown 
457af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4589566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
459af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4609566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4619566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4629566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4639566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
464af33a6ddSJed Brown   }
4659566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
466af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
467af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
468af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
469af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
4709566063dSJacob Faibussowitsch       if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4719566063dSJacob Faibussowitsch       else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
472bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
473af33a6ddSJed Brown     }
474af33a6ddSJed Brown   }
4759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
476af33a6ddSJed Brown 
4779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4789566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
480af33a6ddSJed Brown 
481af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
4829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
48363a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
484af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
485af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
486af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
487af33a6ddSJed Brown 
488af33a6ddSJed Brown     /* for debugging purposes */
489af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
4909371c9d4SSatish Balay       PetscScalar im = interpIndices[0];
4919371c9d4SSatish Balay       PetscScalar jm = interpIndices[1];
492af33a6ddSJed Brown 
49363a3b9bcSJacob 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));
494af33a6ddSJed Brown     }
495af33a6ddSJed Brown 
4969566063dSJacob Faibussowitsch     if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4979566063dSJacob Faibussowitsch     else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
498bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
499af33a6ddSJed Brown   }
5009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
501af33a6ddSJed Brown 
5029566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
5039566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
5049566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
505af33a6ddSJed Brown   if (c->fieldInterpLocal) {
5069566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
5079566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
508af33a6ddSJed Brown   }
5099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
510af33a6ddSJed Brown 
511af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
5129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5139566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5149566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
515af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
516af33a6ddSJed Brown     Qi = c->queue[n];
517bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp];
518af33a6ddSJed Brown   }
5199566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5219566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
522af33a6ddSJed Brown 
523af33a6ddSJed Brown   /* Cleanup */
5249566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpIndices));
5259566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValues));
5269566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValuesOld));
5279566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldValues));
528af33a6ddSJed Brown   PetscFunctionReturn(0);
529af33a6ddSJed Brown }
530af33a6ddSJed Brown 
531*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
532*d71ae5a4SJacob Faibussowitsch {
533af33a6ddSJed Brown   PetscFunctionBegin;
534af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5359566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->neighbors));
5369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5379566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
538af33a6ddSJed Brown   PetscFunctionReturn(0);
539af33a6ddSJed Brown }
540af33a6ddSJed Brown 
541*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
542*d71ae5a4SJacob Faibussowitsch {
543af33a6ddSJed Brown   PetscFunctionBegin;
54463a3b9bcSJacob Faibussowitsch   PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
545af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
546af33a6ddSJed Brown   PetscFunctionReturn(0);
547af33a6ddSJed Brown }
548af33a6ddSJed Brown 
549*d71ae5a4SJacob Faibussowitsch int CharacteristicSendCoordinatesBegin(Characteristic c)
550*d71ae5a4SJacob Faibussowitsch {
551af33a6ddSJed Brown   PetscMPIInt rank, tag = 121;
552af33a6ddSJed Brown   PetscInt    i, n;
553af33a6ddSJed Brown 
554af33a6ddSJed Brown   PetscFunctionBegin;
5559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5569566063dSJacob Faibussowitsch   PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5579566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
558bbd56ea5SKarl Rupp   for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
559af33a6ddSJed Brown   c->fillCount[0] = 0;
56048a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Irecv(&(c->fillCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n - 1])));
56148a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
5629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
563af33a6ddSJed Brown   /* Initialize the remote queue */
564af33a6ddSJed Brown   c->queueLocalMax = c->localOffsets[0] = 0;
565af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
566af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
567af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
568af33a6ddSJed Brown     c->queueRemoteMax += c->fillCount[n];
569af33a6ddSJed Brown     c->localOffsets[n] = c->queueLocalMax;
570af33a6ddSJed Brown     c->queueLocalMax += c->needCount[n];
571af33a6ddSJed Brown   }
572af33a6ddSJed Brown   /* HACK BEGIN */
573bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
574af33a6ddSJed Brown   c->needCount[0] = 0;
575af33a6ddSJed Brown   /* HACK END */
576af33a6ddSJed Brown   if (c->queueRemoteMax) {
5779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5780298fd71SBarry Smith   } else c->queueRemote = NULL;
579af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
580af33a6ddSJed Brown 
581af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
582af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58363a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
5849566063dSJacob 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])));
585af33a6ddSJed Brown   }
586af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
5889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
589af33a6ddSJed Brown   }
590af33a6ddSJed Brown   PetscFunctionReturn(0);
591af33a6ddSJed Brown }
592af33a6ddSJed Brown 
593*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
594*d71ae5a4SJacob Faibussowitsch {
595af33a6ddSJed Brown #if 0
596af33a6ddSJed Brown   PetscMPIInt rank;
597af33a6ddSJed Brown   PetscInt    n;
598af33a6ddSJed Brown #endif
599af33a6ddSJed Brown 
600af33a6ddSJed Brown   PetscFunctionBegin;
6019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
602af33a6ddSJed Brown #if 0
6039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
604af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
6053c633725SBarry 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);
606af33a6ddSJed Brown   }
607af33a6ddSJed Brown #endif
608af33a6ddSJed Brown   PetscFunctionReturn(0);
609af33a6ddSJed Brown }
610af33a6ddSJed Brown 
611*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
612*d71ae5a4SJacob Faibussowitsch {
613af33a6ddSJed Brown   PetscMPIInt tag = 121;
614af33a6ddSJed Brown   PetscInt    n;
615af33a6ddSJed Brown 
616af33a6ddSJed Brown   PetscFunctionBegin;
617af33a6ddSJed Brown   /* SEND AND RECIEVE FILLED REQUESTS for velocities at t_n+1/2 */
61848a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Irecv(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c), &(c->request[n - 1])));
61948a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->queueRemote[c->remoteOffsets[n]]), c->fillCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
620af33a6ddSJed Brown   PetscFunctionReturn(0);
621af33a6ddSJed Brown }
622af33a6ddSJed Brown 
623*d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
624*d71ae5a4SJacob Faibussowitsch {
625af33a6ddSJed Brown   PetscFunctionBegin;
6269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
627af33a6ddSJed Brown   /* Free queue of requests from other procs */
6289566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->queueRemote));
629af33a6ddSJed Brown   PetscFunctionReturn(0);
630af33a6ddSJed Brown }
631af33a6ddSJed Brown 
632af33a6ddSJed Brown /*---------------------------------------------------------------------*/
633af33a6ddSJed Brown /*
634af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
635af33a6ddSJed Brown */
6360b8f38c8SBarry Smith PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
637af33a6ddSJed Brown /*---------------------------------------------------------------------*/
638af33a6ddSJed Brown {
639af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6400b8f38c8SBarry Smith   PetscInt                n;
641af33a6ddSJed Brown 
6420b8f38c8SBarry Smith   PetscFunctionBegin;
643af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6449566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
64548a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
646af33a6ddSJed Brown   }
647af33a6ddSJed Brown 
648af33a6ddSJed Brown   /* SORTING PHASE */
6499371c9d4SSatish Balay   for (n = (size / 2) - 1; n >= 0; n--) { PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/ }
650af33a6ddSJed Brown   for (n = size - 1; n >= 1; n--) {
651af33a6ddSJed Brown     temp     = queue[0];
652af33a6ddSJed Brown     queue[0] = queue[n];
653af33a6ddSJed Brown     queue[n] = temp;
6549566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1));
655af33a6ddSJed Brown   }
656af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
6579566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Avter  Heap sort\n"));
65848a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
659af33a6ddSJed Brown   }
6600b8f38c8SBarry Smith   PetscFunctionReturn(0);
661af33a6ddSJed Brown }
662af33a6ddSJed Brown 
663af33a6ddSJed Brown /*---------------------------------------------------------------------*/
664af33a6ddSJed Brown /*
665af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
666af33a6ddSJed Brown */
6670b8f38c8SBarry Smith PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
668af33a6ddSJed Brown /*---------------------------------------------------------------------*/
669af33a6ddSJed Brown {
670af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
671af33a6ddSJed Brown   PetscInt                maxChild;
672af33a6ddSJed Brown   CharacteristicPointDA2D temp;
673af33a6ddSJed Brown 
674af33a6ddSJed Brown   PetscFunctionBegin;
675af33a6ddSJed Brown   while ((root * 2 <= bottom) && (!done)) {
676af33a6ddSJed Brown     if (root * 2 == bottom) maxChild = root * 2;
677af33a6ddSJed Brown     else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2;
678af33a6ddSJed Brown     else maxChild = root * 2 + 1;
679af33a6ddSJed Brown 
680af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
681af33a6ddSJed Brown       temp            = queue[root];
682af33a6ddSJed Brown       queue[root]     = queue[maxChild];
683af33a6ddSJed Brown       queue[maxChild] = temp;
684af33a6ddSJed Brown       root            = maxChild;
685db4deed7SKarl Rupp     } else done = PETSC_TRUE;
686af33a6ddSJed Brown   }
687af33a6ddSJed Brown   PetscFunctionReturn(0);
688af33a6ddSJed Brown }
689af33a6ddSJed Brown 
690af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
691*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
692*d71ae5a4SJacob Faibussowitsch {
693bff4a2f0SMatthew G. Knepley   DMBoundaryType bx, by;
694af33a6ddSJed Brown   PetscBool      IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
695af33a6ddSJed Brown   MPI_Comm       comm;
696af33a6ddSJed Brown   PetscMPIInt    rank;
697af33a6ddSJed Brown   PetscInt     **procs, pi, pj, pim, pip, pjm, pjp, PI, PJ;
698af33a6ddSJed Brown 
699af33a6ddSJed Brown   PetscFunctionBegin;
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
7019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL));
703af33a6ddSJed Brown 
704bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
705bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
706af33a6ddSJed Brown 
707af33a6ddSJed Brown   neighbors[0] = rank;
708af33a6ddSJed Brown   rank         = 0;
7099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PJ, &procs));
710af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
7119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PI, &(procs[pj])));
712af33a6ddSJed Brown     for (pi = 0; pi < PI; pi++) {
713af33a6ddSJed Brown       procs[pj][pi] = rank;
714af33a6ddSJed Brown       rank++;
715af33a6ddSJed Brown     }
716af33a6ddSJed Brown   }
717af33a6ddSJed Brown 
718af33a6ddSJed Brown   pi  = neighbors[0] % PI;
719af33a6ddSJed Brown   pj  = neighbors[0] / PI;
7209371c9d4SSatish Balay   pim = pi - 1;
7219371c9d4SSatish Balay   if (pim < 0) pim = PI - 1;
722af33a6ddSJed Brown   pip = (pi + 1) % PI;
7239371c9d4SSatish Balay   pjm = pj - 1;
7249371c9d4SSatish Balay   if (pjm < 0) pjm = PJ - 1;
725af33a6ddSJed Brown   pjp = (pj + 1) % PJ;
726af33a6ddSJed Brown 
727af33a6ddSJed Brown   neighbors[1] = procs[pj][pim];
728af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
729af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
730af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
731af33a6ddSJed Brown   neighbors[5] = procs[pj][pip];
732af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
733af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
734af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
735af33a6ddSJed Brown 
736af33a6ddSJed Brown   if (!IPeriodic) {
737af33a6ddSJed Brown     if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0];
738af33a6ddSJed Brown     if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0];
739af33a6ddSJed Brown   }
740af33a6ddSJed Brown 
741af33a6ddSJed Brown   if (!JPeriodic) {
742af33a6ddSJed Brown     if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0];
743af33a6ddSJed Brown     if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0];
744af33a6ddSJed Brown   }
745af33a6ddSJed Brown 
74648a46eb9SPierre Jolivet   for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj]));
7479566063dSJacob Faibussowitsch   PetscCall(PetscFree(procs));
748af33a6ddSJed Brown   PetscFunctionReturn(0);
749af33a6ddSJed Brown }
750af33a6ddSJed Brown 
751af33a6ddSJed Brown /*
752af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
753af33a6ddSJed Brown     2 | 3 | 4
754af33a6ddSJed Brown     __|___|__
755af33a6ddSJed Brown     1 | 0 | 5
756af33a6ddSJed Brown     __|___|__
757af33a6ddSJed Brown     8 | 7 | 6
758af33a6ddSJed Brown       |   |
759af33a6ddSJed Brown */
760*d71ae5a4SJacob Faibussowitsch PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
761*d71ae5a4SJacob Faibussowitsch {
762af33a6ddSJed Brown   DMDALocalInfo info;
7631cc8b206SBarry Smith   PetscReal     is, ie, js, je;
764af33a6ddSJed Brown 
7659566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
7669371c9d4SSatish Balay   is = (PetscReal)info.xs - 0.5;
7679371c9d4SSatish Balay   ie = (PetscReal)info.xs + info.xm - 0.5;
7689371c9d4SSatish Balay   js = (PetscReal)info.ys - 0.5;
7699371c9d4SSatish Balay   je = (PetscReal)info.ys + info.ym - 0.5;
770af33a6ddSJed Brown 
771af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
772bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
773bbd56ea5SKarl Rupp     else if (jr < js) return 7;
774bbd56ea5SKarl Rupp     else return 3;
775af33a6ddSJed Brown   } else if (ir < is) { /* left column */
776bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
777bbd56ea5SKarl Rupp     else if (jr < js) return 8;
778bbd56ea5SKarl Rupp     else return 2;
779af33a6ddSJed Brown   } else { /* right column */
780bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
781bbd56ea5SKarl Rupp     else if (jr < js) return 6;
782bbd56ea5SKarl Rupp     else return 4;
783af33a6ddSJed Brown   }
784af33a6ddSJed Brown }
785