xref: /petsc/src/ts/characteristic/interface/characteristic.c (revision 66976f2f44dcc61d86a452a70219fb23b45d00f0)
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 
16*66976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM, PetscMPIInt[]);
17*66976f2fSJacob Faibussowitsch static PetscInt       DMDAGetNeighborRelative(DM, PetscReal, PetscReal);
18af33a6ddSJed Brown 
19*66976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic, Queue, PetscInt);
20*66976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic, Queue, PetscInt, PetscInt);
21af33a6ddSJed Brown 
22*66976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicView(Characteristic c, PetscViewer viewer)
23d71ae5a4SJacob Faibussowitsch {
24af33a6ddSJed Brown   PetscBool iascii;
25af33a6ddSJed Brown 
26af33a6ddSJed Brown   PetscFunctionBegin;
27af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
2848a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)c), &viewer));
29af33a6ddSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
30af33a6ddSJed Brown   PetscCheckSameComm(c, 1, viewer, 2);
31af33a6ddSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
33ad540459SPierre Jolivet   if (!iascii) PetscTryTypeMethod(c, view, viewer);
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35af33a6ddSJed Brown }
36af33a6ddSJed Brown 
37d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicDestroy(Characteristic *c)
38d71ae5a4SJacob Faibussowitsch {
39af33a6ddSJed Brown   PetscFunctionBegin;
403ba16761SJacob Faibussowitsch   if (!*c) PetscFunctionReturn(PETSC_SUCCESS);
416bf464f9SBarry Smith   PetscValidHeaderSpecific(*c, CHARACTERISTIC_CLASSID, 1);
423ba16761SJacob Faibussowitsch   if (--((PetscObject)(*c))->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
43af33a6ddSJed Brown 
4448a46eb9SPierre Jolivet   if ((*c)->ops->destroy) PetscCall((*(*c)->ops->destroy)((*c)));
459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&(*c)->itemType));
469566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queue));
479566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueLocal));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->queueRemote));
499566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->neighbors));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->needCount));
519566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->localOffsets));
529566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->fillCount));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->remoteOffsets));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->request));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree((*c)->status));
569566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(c));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58af33a6ddSJed Brown }
59af33a6ddSJed Brown 
60d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicCreate(MPI_Comm comm, Characteristic *c)
61d71ae5a4SJacob Faibussowitsch {
62af33a6ddSJed Brown   Characteristic newC;
63af33a6ddSJed Brown 
64af33a6ddSJed Brown   PetscFunctionBegin;
654f572ea9SToby Isaac   PetscAssertPointer(c, 2);
660298fd71SBarry Smith   *c = NULL;
679566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
68af33a6ddSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(newC, CHARACTERISTIC_CLASSID, "Characteristic", "Characteristic", "Characteristic", comm, CharacteristicDestroy, CharacteristicView));
70af33a6ddSJed Brown   *c = newC;
71af33a6ddSJed Brown 
72af33a6ddSJed Brown   newC->structured          = PETSC_TRUE;
73af33a6ddSJed Brown   newC->numIds              = 0;
740298fd71SBarry Smith   newC->velocityDA          = NULL;
750298fd71SBarry Smith   newC->velocity            = NULL;
760298fd71SBarry Smith   newC->velocityOld         = NULL;
77af33a6ddSJed Brown   newC->numVelocityComp     = 0;
780298fd71SBarry Smith   newC->velocityComp        = NULL;
790298fd71SBarry Smith   newC->velocityInterp      = NULL;
800298fd71SBarry Smith   newC->velocityInterpLocal = NULL;
810298fd71SBarry Smith   newC->velocityCtx         = NULL;
820298fd71SBarry Smith   newC->fieldDA             = NULL;
830298fd71SBarry Smith   newC->field               = NULL;
84af33a6ddSJed Brown   newC->numFieldComp        = 0;
850298fd71SBarry Smith   newC->fieldComp           = NULL;
860298fd71SBarry Smith   newC->fieldInterp         = NULL;
870298fd71SBarry Smith   newC->fieldInterpLocal    = NULL;
880298fd71SBarry Smith   newC->fieldCtx            = NULL;
890298fd71SBarry Smith   newC->itemType            = 0;
900298fd71SBarry Smith   newC->queue               = NULL;
91af33a6ddSJed Brown   newC->queueSize           = 0;
92af33a6ddSJed Brown   newC->queueMax            = 0;
930298fd71SBarry Smith   newC->queueLocal          = NULL;
94af33a6ddSJed Brown   newC->queueLocalSize      = 0;
95af33a6ddSJed Brown   newC->queueLocalMax       = 0;
960298fd71SBarry Smith   newC->queueRemote         = NULL;
97af33a6ddSJed Brown   newC->queueRemoteSize     = 0;
98af33a6ddSJed Brown   newC->queueRemoteMax      = 0;
99af33a6ddSJed Brown   newC->numNeighbors        = 0;
1000298fd71SBarry Smith   newC->neighbors           = NULL;
1010298fd71SBarry Smith   newC->needCount           = NULL;
1020298fd71SBarry Smith   newC->localOffsets        = NULL;
1030298fd71SBarry Smith   newC->fillCount           = NULL;
1040298fd71SBarry Smith   newC->remoteOffsets       = NULL;
1050298fd71SBarry Smith   newC->request             = NULL;
1060298fd71SBarry Smith   newC->status              = NULL;
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108af33a6ddSJed Brown }
109af33a6ddSJed Brown 
110af33a6ddSJed Brown /*@C
111af33a6ddSJed Brown   CharacteristicSetType - Builds Characteristic for a particular solver.
112af33a6ddSJed Brown 
113c3339decSBarry Smith   Logically Collective
114af33a6ddSJed Brown 
115af33a6ddSJed Brown   Input Parameters:
116af33a6ddSJed Brown + c    - the method of characteristics context
117af33a6ddSJed Brown - type - a known method
118af33a6ddSJed Brown 
119af33a6ddSJed Brown   Options Database Key:
120af33a6ddSJed Brown . -characteristic_type <method> - Sets the method; use -help for a list
121af33a6ddSJed Brown     of available methods
122af33a6ddSJed Brown 
123bcf0153eSBarry Smith   Level: intermediate
124bcf0153eSBarry Smith 
125af33a6ddSJed Brown   Notes:
12630a76a96SBarry Smith   See "include/petsccharacteristic.h" for available methods
127af33a6ddSJed Brown 
128af33a6ddSJed Brown   Normally, it is best to use the CharacteristicSetFromOptions() command and
129af33a6ddSJed Brown   then set the Characteristic type from the options database rather than by using
130af33a6ddSJed Brown   this routine.  Using the options database provides the user with
131af33a6ddSJed Brown   maximum flexibility in evaluating the many different Krylov methods.
132af33a6ddSJed Brown   The CharacteristicSetType() routine is provided for those situations where it
133af33a6ddSJed Brown   is necessary to set the iterative solver independently of the command
134af33a6ddSJed Brown   line or options database.  This might be the case, for example, when
135af33a6ddSJed Brown   the choice of iterative solver changes during the execution of the
136af33a6ddSJed Brown   program, and the user's application is taking responsibility for
137af33a6ddSJed Brown   choosing the appropriate method.  In other words, this routine is
138af33a6ddSJed Brown   not for beginners.
139af33a6ddSJed Brown 
1401cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicType`
141af33a6ddSJed Brown @*/
142d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetType(Characteristic c, CharacteristicType type)
143d71ae5a4SJacob Faibussowitsch {
144af33a6ddSJed Brown   PetscBool match;
1455f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(Characteristic);
146af33a6ddSJed Brown 
147af33a6ddSJed Brown   PetscFunctionBegin;
148af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
1494f572ea9SToby Isaac   PetscAssertPointer(type, 2);
150af33a6ddSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)c, type, &match));
1523ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
153af33a6ddSJed Brown 
154af33a6ddSJed Brown   if (c->data) {
155af33a6ddSJed Brown     /* destroy the old private Characteristic context */
156dbbe0bcdSBarry Smith     PetscUseTypeMethod(c, destroy);
1570298fd71SBarry Smith     c->ops->destroy = NULL;
1582120983fSLisandro Dalcin     c->data         = NULL;
159af33a6ddSJed Brown   }
160af33a6ddSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(CharacteristicList, type, &r));
1626adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)c), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown Characteristic type given: %s", type);
163af33a6ddSJed Brown   c->setupcalled = 0;
1649566063dSJacob Faibussowitsch   PetscCall((*r)(c));
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)c, type));
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167af33a6ddSJed Brown }
168af33a6ddSJed Brown 
169af33a6ddSJed Brown /*@
170af33a6ddSJed Brown   CharacteristicSetUp - Sets up the internal data structures for the
171af33a6ddSJed Brown   later use of an iterative solver.
172af33a6ddSJed Brown 
173c3339decSBarry Smith   Collective
174af33a6ddSJed Brown 
175af33a6ddSJed Brown   Input Parameter:
176b43aa488SJacob Faibussowitsch . c - iterative context obtained from CharacteristicCreate()
177af33a6ddSJed Brown 
178af33a6ddSJed Brown   Level: developer
179af33a6ddSJed Brown 
1801cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicCreate()`, `CharacteristicSolve()`, `CharacteristicDestroy()`
181af33a6ddSJed Brown @*/
182d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetUp(Characteristic c)
183d71ae5a4SJacob Faibussowitsch {
184af33a6ddSJed Brown   PetscFunctionBegin;
185af33a6ddSJed Brown   PetscValidHeaderSpecific(c, CHARACTERISTIC_CLASSID, 1);
186af33a6ddSJed Brown 
18748a46eb9SPierre Jolivet   if (!((PetscObject)c)->type_name) PetscCall(CharacteristicSetType(c, CHARACTERISTICDA));
188af33a6ddSJed Brown 
1893ba16761SJacob Faibussowitsch   if (c->setupcalled == 2) PetscFunctionReturn(PETSC_SUCCESS);
190af33a6ddSJed Brown 
1919566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
192ad540459SPierre Jolivet   if (!c->setupcalled) PetscUseTypeMethod(c, setup);
1939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_SetUp, c, NULL, NULL, NULL));
194af33a6ddSJed Brown   c->setupcalled = 2;
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
196af33a6ddSJed Brown }
197af33a6ddSJed Brown 
198af33a6ddSJed Brown /*@C
1991c84c290SBarry Smith   CharacteristicRegister -  Adds a solver to the method of characteristics package.
2001c84c290SBarry Smith 
2011c84c290SBarry Smith   Not Collective
2021c84c290SBarry Smith 
2031c84c290SBarry Smith   Input Parameters:
2042fe279fdSBarry Smith + sname    - name of a new user-defined solver
2052fe279fdSBarry Smith - function - routine to create method context
2061c84c290SBarry Smith 
207bcf0153eSBarry Smith   Level: advanced
208bcf0153eSBarry Smith 
209b43aa488SJacob Faibussowitsch   Example 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 
2271cc06b55SBarry Smith .seealso: [](ch_ts), `CharacteristicRegisterAll()`, `CharacteristicRegisterDestroy()`
228af33a6ddSJed Brown @*/
229d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicRegister(const char sname[], PetscErrorCode (*function)(Characteristic))
230d71ae5a4SJacob Faibussowitsch {
231af33a6ddSJed Brown   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(CharacteristicInitializePackage());
2339566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&CharacteristicList, sname, function));
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235af33a6ddSJed Brown }
236af33a6ddSJed Brown 
237d71ae5a4SJacob 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)
238d71ae5a4SJacob Faibussowitsch {
239af33a6ddSJed Brown   PetscFunctionBegin;
240af33a6ddSJed Brown   c->velocityDA      = da;
241af33a6ddSJed Brown   c->velocity        = v;
242af33a6ddSJed Brown   c->velocityOld     = vOld;
243af33a6ddSJed Brown   c->numVelocityComp = numComponents;
244af33a6ddSJed Brown   c->velocityComp    = components;
245af33a6ddSJed Brown   c->velocityInterp  = interp;
246af33a6ddSJed Brown   c->velocityCtx     = ctx;
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
248af33a6ddSJed Brown }
249af33a6ddSJed Brown 
250d71ae5a4SJacob 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)
251d71ae5a4SJacob Faibussowitsch {
252af33a6ddSJed Brown   PetscFunctionBegin;
253af33a6ddSJed Brown   c->velocityDA          = da;
254af33a6ddSJed Brown   c->velocity            = v;
255af33a6ddSJed Brown   c->velocityOld         = vOld;
256af33a6ddSJed Brown   c->numVelocityComp     = numComponents;
257af33a6ddSJed Brown   c->velocityComp        = components;
258af33a6ddSJed Brown   c->velocityInterpLocal = interp;
259af33a6ddSJed Brown   c->velocityCtx         = ctx;
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
261af33a6ddSJed Brown }
262af33a6ddSJed Brown 
263d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolation(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(Vec, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
264d71ae5a4SJacob Faibussowitsch {
265af33a6ddSJed Brown   PetscFunctionBegin;
266af33a6ddSJed Brown #if 0
2673c633725SBarry 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.");
268af33a6ddSJed Brown #endif
269af33a6ddSJed Brown   c->fieldDA      = da;
270af33a6ddSJed Brown   c->field        = v;
271af33a6ddSJed Brown   c->numFieldComp = numComponents;
272af33a6ddSJed Brown   c->fieldComp    = components;
273af33a6ddSJed Brown   c->fieldInterp  = interp;
274af33a6ddSJed Brown   c->fieldCtx     = ctx;
2753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276af33a6ddSJed Brown }
277af33a6ddSJed Brown 
278d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetFieldInterpolationLocal(Characteristic c, DM da, Vec v, PetscInt numComponents, PetscInt components[], PetscErrorCode (*interp)(void *, PetscReal[], PetscInt, PetscInt[], PetscScalar[], void *), void *ctx)
279d71ae5a4SJacob Faibussowitsch {
280af33a6ddSJed Brown   PetscFunctionBegin;
281af33a6ddSJed Brown #if 0
2823c633725SBarry 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.");
283af33a6ddSJed Brown #endif
284af33a6ddSJed Brown   c->fieldDA          = da;
285af33a6ddSJed Brown   c->field            = v;
286af33a6ddSJed Brown   c->numFieldComp     = numComponents;
287af33a6ddSJed Brown   c->fieldComp        = components;
288af33a6ddSJed Brown   c->fieldInterpLocal = interp;
289af33a6ddSJed Brown   c->fieldCtx         = ctx;
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
291af33a6ddSJed Brown }
292af33a6ddSJed Brown 
293d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSolve(Characteristic c, PetscReal dt, Vec solution)
294d71ae5a4SJacob Faibussowitsch {
295af33a6ddSJed Brown   CharacteristicPointDA2D Qi;
296af33a6ddSJed Brown   DM                      da = c->velocityDA;
297af33a6ddSJed Brown   Vec                     velocityLocal, velocityLocalOld;
298af33a6ddSJed Brown   Vec                     fieldLocal;
299af33a6ddSJed Brown   DMDALocalInfo           info;
300af33a6ddSJed Brown   PetscScalar           **solArray;
301af33a6ddSJed Brown   void                   *velocityArray;
302af33a6ddSJed Brown   void                   *velocityArrayOld;
303af33a6ddSJed Brown   void                   *fieldArray;
3041cc8b206SBarry Smith   PetscScalar            *interpIndices;
3051cc8b206SBarry Smith   PetscScalar            *velocityValues, *velocityValuesOld;
3061cc8b206SBarry Smith   PetscScalar            *fieldValues;
307af33a6ddSJed Brown   PetscMPIInt             rank;
308af33a6ddSJed Brown   PetscInt                dim;
309af33a6ddSJed Brown   PetscMPIInt             neighbors[9];
310af33a6ddSJed Brown   PetscInt                dof;
311af33a6ddSJed Brown   PetscInt                gx, gy;
312af33a6ddSJed Brown   PetscInt                n, is, ie, js, je, comp;
313af33a6ddSJed Brown 
314af33a6ddSJed Brown   PetscFunctionBegin;
315af33a6ddSJed Brown   c->queueSize = 0;
3169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
3179566063dSJacob Faibussowitsch   PetscCall(DMDAGetNeighborsRank(da, neighbors));
3189566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetNeighbors(c, 9, neighbors));
3199566063dSJacob Faibussowitsch   PetscCall(CharacteristicSetUp(c));
320af33a6ddSJed Brown   /* global and local grid info */
3219566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &gx, &gy, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
3239371c9d4SSatish Balay   is = info.xs;
3249371c9d4SSatish Balay   ie = info.xs + info.xm;
3259371c9d4SSatish Balay   js = info.ys;
3269371c9d4SSatish Balay   je = info.ys + info.ym;
327af33a6ddSJed Brown   /* Allocation */
3289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim, &interpIndices));
3299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValues));
3309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numVelocityComp, &velocityValuesOld));
3319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->numFieldComp, &fieldValues));
3329566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
333af33a6ddSJed Brown 
334af33a6ddSJed Brown   /* -----------------------------------------------------------------------
335af33a6ddSJed Brown      PART 1, AT t-dt/2
336af33a6ddSJed Brown      -----------------------------------------------------------------------*/
3379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
338af33a6ddSJed Brown   /* GET POSITION AT HALF TIME IN THE PAST */
339af33a6ddSJed Brown   if (c->velocityInterpLocal) {
3409566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocal));
3419566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->velocityDA, &velocityLocalOld));
3429566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3439566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocity, INSERT_VALUES, velocityLocal));
3449566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3459566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->velocityDA, c->velocityOld, INSERT_VALUES, velocityLocalOld));
3469566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocal, &velocityArray));
3479566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
348af33a6ddSJed Brown   }
3499566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n - 1/2}\n"));
350af33a6ddSJed Brown   for (Qi.j = js; Qi.j < je; Qi.j++) {
351af33a6ddSJed Brown     for (Qi.i = is; Qi.i < ie; Qi.i++) {
352af33a6ddSJed Brown       interpIndices[0] = Qi.i;
353af33a6ddSJed Brown       interpIndices[1] = Qi.j;
3549566063dSJacob Faibussowitsch       if (c->velocityInterpLocal) PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3559566063dSJacob Faibussowitsch       else PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
356af33a6ddSJed Brown       Qi.x = Qi.i - velocityValues[0] * dt / 2.0;
357af33a6ddSJed Brown       Qi.y = Qi.j - velocityValues[1] * dt / 2.0;
358af33a6ddSJed Brown 
359af33a6ddSJed Brown       /* Determine whether the position at t - dt/2 is local */
360af33a6ddSJed Brown       Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
361af33a6ddSJed Brown 
362af33a6ddSJed Brown       /* Check for Periodic boundaries and move all periodic points back onto the domain */
3639566063dSJacob Faibussowitsch       PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
3649566063dSJacob Faibussowitsch       PetscCall(CharacteristicAddPoint(c, &Qi));
365af33a6ddSJed Brown     }
366af33a6ddSJed Brown   }
3679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_QueueSetup, NULL, NULL, NULL, NULL));
368af33a6ddSJed Brown 
3699566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3709566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
3719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
372af33a6ddSJed Brown 
3739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
374af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (local values) */
3759566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local velocities at t_{n - 1/2}\n"));
376af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
377af33a6ddSJed Brown     Qi = c->queue[n];
378af33a6ddSJed Brown     if (c->neighbors[Qi.proc] == rank) {
379af33a6ddSJed Brown       interpIndices[0] = Qi.x;
380af33a6ddSJed Brown       interpIndices[1] = Qi.y;
381af33a6ddSJed Brown       if (c->velocityInterpLocal) {
3829566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3839566063dSJacob Faibussowitsch         PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
384af33a6ddSJed Brown       } else {
3859566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
3869566063dSJacob Faibussowitsch         PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
387af33a6ddSJed Brown       }
388af33a6ddSJed Brown       Qi.x = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
389af33a6ddSJed Brown       Qi.y = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
390af33a6ddSJed Brown     }
391af33a6ddSJed Brown     c->queue[n] = Qi;
392af33a6ddSJed Brown   }
3939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeLocal, NULL, NULL, NULL, NULL));
394af33a6ddSJed Brown 
3959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
3969566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
3979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
398af33a6ddSJed Brown 
399af33a6ddSJed Brown   /* Calculate velocity at t_n+1/2 (fill remote requests) */
4009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
40163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote velocities at t_{n - 1/2}\n", c->queueRemoteSize));
402af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
403af33a6ddSJed Brown     Qi               = c->queueRemote[n];
404af33a6ddSJed Brown     interpIndices[0] = Qi.x;
405af33a6ddSJed Brown     interpIndices[1] = Qi.y;
406af33a6ddSJed Brown     if (c->velocityInterpLocal) {
4079566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArray, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4089566063dSJacob Faibussowitsch       PetscCall(c->velocityInterpLocal(velocityArrayOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
409af33a6ddSJed Brown     } else {
4109566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocity, interpIndices, c->numVelocityComp, c->velocityComp, velocityValues, c->velocityCtx));
4119566063dSJacob Faibussowitsch       PetscCall(c->velocityInterp(c->velocityOld, interpIndices, c->numVelocityComp, c->velocityComp, velocityValuesOld, c->velocityCtx));
412af33a6ddSJed Brown     }
413af33a6ddSJed Brown     Qi.x              = 0.5 * (velocityValues[0] + velocityValuesOld[0]);
414af33a6ddSJed Brown     Qi.y              = 0.5 * (velocityValues[1] + velocityValuesOld[1]);
415af33a6ddSJed Brown     c->queueRemote[n] = Qi;
416af33a6ddSJed Brown   }
4179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeRemote, NULL, NULL, NULL, NULL));
4189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
4199566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
4209566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
421af33a6ddSJed Brown   if (c->velocityInterpLocal) {
4229566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocal, &velocityArray));
4239566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->velocityDA, velocityLocalOld, &velocityArrayOld));
4249566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocal));
4259566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->velocityDA, &velocityLocalOld));
426af33a6ddSJed Brown   }
4279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_HalfTimeExchange, NULL, NULL, NULL, NULL));
428af33a6ddSJed Brown 
429af33a6ddSJed Brown   /* -----------------------------------------------------------------------
430af33a6ddSJed Brown      PART 2, AT t-dt
431af33a6ddSJed Brown      -----------------------------------------------------------------------*/
432af33a6ddSJed Brown 
433af33a6ddSJed Brown   /* GET POSITION AT t_n (local values) */
4349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
4359566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating position at t_{n}\n"));
436af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
437af33a6ddSJed Brown     Qi   = c->queue[n];
438af33a6ddSJed Brown     Qi.x = Qi.i - Qi.x * dt;
439af33a6ddSJed Brown     Qi.y = Qi.j - Qi.y * dt;
440af33a6ddSJed Brown 
441af33a6ddSJed Brown     /* Determine whether the position at t-dt is local */
442af33a6ddSJed Brown     Qi.proc = DMDAGetNeighborRelative(da, Qi.x, Qi.y);
443af33a6ddSJed Brown 
444af33a6ddSJed Brown     /* Check for Periodic boundaries and move all periodic points back onto the domain */
4459566063dSJacob Faibussowitsch     PetscCall(DMDAMapCoordsToPeriodicDomain(da, &(Qi.x), &(Qi.y)));
446af33a6ddSJed Brown 
447af33a6ddSJed Brown     c->queue[n] = Qi;
448af33a6ddSJed Brown   }
4499566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
450af33a6ddSJed Brown 
4519566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4529566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesBegin(c));
4539566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
454af33a6ddSJed Brown 
455af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (LOCAL REQUESTS) */
4569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
457af33a6ddSJed Brown   if (c->fieldInterpLocal) {
4589566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(c->fieldDA, &fieldLocal));
4599566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4609566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(c->fieldDA, c->field, INSERT_VALUES, fieldLocal));
4619566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(c->fieldDA, fieldLocal, &fieldArray));
462af33a6ddSJed Brown   }
4639566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating local field at t_{n}\n"));
464af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
465af33a6ddSJed Brown     if (c->neighbors[c->queue[n].proc] == rank) {
466af33a6ddSJed Brown       interpIndices[0] = c->queue[n].x;
467af33a6ddSJed Brown       interpIndices[1] = c->queue[n].y;
4689566063dSJacob Faibussowitsch       if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4699566063dSJacob Faibussowitsch       else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
470bbd56ea5SKarl Rupp       for (comp = 0; comp < c->numFieldComp; comp++) c->queue[n].field[comp] = fieldValues[comp];
471af33a6ddSJed Brown     }
472af33a6ddSJed Brown   }
4739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeLocal, NULL, NULL, NULL, NULL));
474af33a6ddSJed Brown 
4759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
4769566063dSJacob Faibussowitsch   PetscCall(CharacteristicSendCoordinatesEnd(c));
4779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
478af33a6ddSJed Brown 
479af33a6ddSJed Brown   /* GET VALUE AT FULL TIME IN THE PAST (REMOTE REQUESTS) */
4809566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
48163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Calculating %" PetscInt_FMT " remote field points at t_{n}\n", c->queueRemoteSize));
482af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
483af33a6ddSJed Brown     interpIndices[0] = c->queueRemote[n].x;
484af33a6ddSJed Brown     interpIndices[1] = c->queueRemote[n].y;
485af33a6ddSJed Brown 
486af33a6ddSJed Brown     /* for debugging purposes */
487af33a6ddSJed Brown     if (1) { /* hacked bounds test...let's do better */
4889371c9d4SSatish Balay       PetscScalar im = interpIndices[0];
4899371c9d4SSatish Balay       PetscScalar jm = interpIndices[1];
490af33a6ddSJed Brown 
49163a3b9bcSJacob 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));
492af33a6ddSJed Brown     }
493af33a6ddSJed Brown 
4949566063dSJacob Faibussowitsch     if (c->fieldInterpLocal) PetscCall(c->fieldInterpLocal(fieldArray, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
4959566063dSJacob Faibussowitsch     else PetscCall(c->fieldInterp(c->field, interpIndices, c->numFieldComp, c->fieldComp, fieldValues, c->fieldCtx));
496bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) c->queueRemote[n].field[comp] = fieldValues[comp];
497af33a6ddSJed Brown   }
4989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeRemote, NULL, NULL, NULL, NULL));
499af33a6ddSJed Brown 
5009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
5019566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesBegin(c));
5029566063dSJacob Faibussowitsch   PetscCall(CharacteristicGetValuesEnd(c));
503af33a6ddSJed Brown   if (c->fieldInterpLocal) {
5049566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(c->fieldDA, fieldLocal, &fieldArray));
5059566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(c->fieldDA, &fieldLocal));
506af33a6ddSJed Brown   }
5079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_FullTimeExchange, NULL, NULL, NULL, NULL));
508af33a6ddSJed Brown 
509af33a6ddSJed Brown   /* Return field of characteristics at t_n-1 */
5109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5119566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(c->fieldDA, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5129566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(c->fieldDA, solution, &solArray));
513af33a6ddSJed Brown   for (n = 0; n < c->queueSize; n++) {
514af33a6ddSJed Brown     Qi = c->queue[n];
515bbd56ea5SKarl Rupp     for (comp = 0; comp < c->numFieldComp; comp++) solArray[Qi.j][Qi.i * dof + c->fieldComp[comp]] = Qi.field[comp];
516af33a6ddSJed Brown   }
5179566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(c->fieldDA, solution, &solArray));
5189566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_DAUpdate, NULL, NULL, NULL, NULL));
5199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(CHARACTERISTIC_Solve, NULL, NULL, NULL, NULL));
520af33a6ddSJed Brown 
521af33a6ddSJed Brown   /* Cleanup */
5229566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpIndices));
5239566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValues));
5249566063dSJacob Faibussowitsch   PetscCall(PetscFree(velocityValuesOld));
5259566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldValues));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527af33a6ddSJed Brown }
528af33a6ddSJed Brown 
529d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSetNeighbors(Characteristic c, PetscInt numNeighbors, PetscMPIInt neighbors[])
530d71ae5a4SJacob Faibussowitsch {
531af33a6ddSJed Brown   PetscFunctionBegin;
532af33a6ddSJed Brown   c->numNeighbors = numNeighbors;
5339566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->neighbors));
5349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &c->neighbors));
5359566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(c->neighbors, neighbors, numNeighbors));
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
537af33a6ddSJed Brown }
538af33a6ddSJed Brown 
539d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicAddPoint(Characteristic c, CharacteristicPointDA2D *point)
540d71ae5a4SJacob Faibussowitsch {
541af33a6ddSJed Brown   PetscFunctionBegin;
54263a3b9bcSJacob Faibussowitsch   PetscCheck(c->queueSize < c->queueMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exceeded maximum queue size %" PetscInt_FMT, c->queueMax);
543af33a6ddSJed Brown   c->queue[c->queueSize++] = *point;
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
545af33a6ddSJed Brown }
546af33a6ddSJed Brown 
5473ba16761SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesBegin(Characteristic c)
548d71ae5a4SJacob Faibussowitsch {
549af33a6ddSJed Brown   PetscMPIInt rank, tag = 121;
550af33a6ddSJed Brown   PetscInt    i, n;
551af33a6ddSJed Brown 
552af33a6ddSJed Brown   PetscFunctionBegin;
5539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
5549566063dSJacob Faibussowitsch   PetscCall(CharacteristicHeapSort(c, c->queue, c->queueSize));
5559566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(c->needCount, c->numNeighbors));
556bbd56ea5SKarl Rupp   for (i = 0; i < c->queueSize; i++) c->needCount[c->queue[i].proc]++;
557af33a6ddSJed Brown   c->fillCount[0] = 0;
55848a46eb9SPierre 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])));
55948a46eb9SPierre Jolivet   for (n = 1; n < c->numNeighbors; n++) PetscCallMPI(MPI_Send(&(c->needCount[n]), 1, MPIU_INT, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
5609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
561af33a6ddSJed Brown   /* Initialize the remote queue */
562af33a6ddSJed Brown   c->queueLocalMax = c->localOffsets[0] = 0;
563af33a6ddSJed Brown   c->queueRemoteMax = c->remoteOffsets[0] = 0;
564af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
565af33a6ddSJed Brown     c->remoteOffsets[n] = c->queueRemoteMax;
566af33a6ddSJed Brown     c->queueRemoteMax += c->fillCount[n];
567af33a6ddSJed Brown     c->localOffsets[n] = c->queueLocalMax;
568af33a6ddSJed Brown     c->queueLocalMax += c->needCount[n];
569af33a6ddSJed Brown   }
570af33a6ddSJed Brown   /* HACK BEGIN */
571bbd56ea5SKarl Rupp   for (n = 1; n < c->numNeighbors; n++) c->localOffsets[n] += c->needCount[0];
572af33a6ddSJed Brown   c->needCount[0] = 0;
573af33a6ddSJed Brown   /* HACK END */
574af33a6ddSJed Brown   if (c->queueRemoteMax) {
5759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(c->queueRemoteMax, &c->queueRemote));
5760298fd71SBarry Smith   } else c->queueRemote = NULL;
577af33a6ddSJed Brown   c->queueRemoteSize = c->queueRemoteMax;
578af33a6ddSJed Brown 
579af33a6ddSJed Brown   /* Send and Receive requests for values at t_n+1/2, giving the coordinates for interpolation */
580af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Receiving %" PetscInt_FMT " requests for values from proc %d\n", c->fillCount[n], c->neighbors[n]));
5829566063dSJacob 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])));
583af33a6ddSJed Brown   }
584af33a6ddSJed Brown   for (n = 1; n < c->numNeighbors; n++) {
58563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Sending %" PetscInt_FMT " requests for values from proc %d\n", c->needCount[n], c->neighbors[n]));
5869566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&(c->queue[c->localOffsets[n]]), c->needCount[n], c->itemType, c->neighbors[n], tag, PetscObjectComm((PetscObject)c)));
587af33a6ddSJed Brown   }
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
589af33a6ddSJed Brown }
590af33a6ddSJed Brown 
591d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicSendCoordinatesEnd(Characteristic c)
592d71ae5a4SJacob Faibussowitsch {
593af33a6ddSJed Brown #if 0
594af33a6ddSJed Brown   PetscMPIInt rank;
595af33a6ddSJed Brown   PetscInt    n;
596af33a6ddSJed Brown #endif
597af33a6ddSJed Brown 
598af33a6ddSJed Brown   PetscFunctionBegin;
5999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
600af33a6ddSJed Brown #if 0
6019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)c), &rank));
602af33a6ddSJed Brown   for (n = 0; n < c->queueRemoteSize; n++) {
6033c633725SBarry 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);
604af33a6ddSJed Brown   }
605af33a6ddSJed Brown #endif
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607af33a6ddSJed Brown }
608af33a6ddSJed Brown 
609d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesBegin(Characteristic c)
610d71ae5a4SJacob Faibussowitsch {
611af33a6ddSJed Brown   PetscMPIInt tag = 121;
612af33a6ddSJed Brown   PetscInt    n;
613af33a6ddSJed Brown 
614af33a6ddSJed Brown   PetscFunctionBegin;
615d5b43468SJose E. Roman   /* SEND AND RECEIVE FILLED REQUESTS for velocities at t_n+1/2 */
61648a46eb9SPierre 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])));
61748a46eb9SPierre 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)));
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619af33a6ddSJed Brown }
620af33a6ddSJed Brown 
621d71ae5a4SJacob Faibussowitsch PetscErrorCode CharacteristicGetValuesEnd(Characteristic c)
622d71ae5a4SJacob Faibussowitsch {
623af33a6ddSJed Brown   PetscFunctionBegin;
6249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(c->numNeighbors - 1, c->request, c->status));
625af33a6ddSJed Brown   /* Free queue of requests from other procs */
6269566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->queueRemote));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
628af33a6ddSJed Brown }
629af33a6ddSJed Brown 
630af33a6ddSJed Brown /*---------------------------------------------------------------------*/
631af33a6ddSJed Brown /*
632af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
633af33a6ddSJed Brown */
634*66976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicHeapSort(Characteristic c, Queue queue, PetscInt size)
635af33a6ddSJed Brown /*---------------------------------------------------------------------*/
636af33a6ddSJed Brown {
637af33a6ddSJed Brown   CharacteristicPointDA2D temp;
6380b8f38c8SBarry Smith   PetscInt                n;
639af33a6ddSJed Brown 
6400b8f38c8SBarry Smith   PetscFunctionBegin;
641af33a6ddSJed Brown   if (0) { /* Check the order of the queue before sorting */
6429566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Before Heap sort\n"));
64348a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
644af33a6ddSJed Brown   }
645af33a6ddSJed Brown 
646af33a6ddSJed Brown   /* SORTING PHASE */
6479371c9d4SSatish Balay   for (n = (size / 2) - 1; n >= 0; n--) { PetscCall(CharacteristicSiftDown(c, queue, n, size - 1)); /* Rich had size-1 here, Matt had size*/ }
648af33a6ddSJed Brown   for (n = size - 1; n >= 1; n--) {
649af33a6ddSJed Brown     temp     = queue[0];
650af33a6ddSJed Brown     queue[0] = queue[n];
651af33a6ddSJed Brown     queue[n] = temp;
6529566063dSJacob Faibussowitsch     PetscCall(CharacteristicSiftDown(c, queue, 0, n - 1));
653af33a6ddSJed Brown   }
654af33a6ddSJed Brown   if (0) { /* Check the order of the queue after sorting */
6559566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Avter  Heap sort\n"));
65648a46eb9SPierre Jolivet     for (n = 0; n < size; n++) PetscCall(PetscInfo(NULL, "%" PetscInt_FMT " %d\n", n, queue[n].proc));
657af33a6ddSJed Brown   }
6583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
659af33a6ddSJed Brown }
660af33a6ddSJed Brown 
661af33a6ddSJed Brown /*---------------------------------------------------------------------*/
662af33a6ddSJed Brown /*
663af33a6ddSJed Brown   Based on code from http://linux.wku.edu/~lamonml/algor/sort/heap.html
664af33a6ddSJed Brown */
665*66976f2fSJacob Faibussowitsch static PetscErrorCode CharacteristicSiftDown(Characteristic c, Queue queue, PetscInt root, PetscInt bottom)
666af33a6ddSJed Brown /*---------------------------------------------------------------------*/
667af33a6ddSJed Brown {
668af33a6ddSJed Brown   PetscBool               done = PETSC_FALSE;
669af33a6ddSJed Brown   PetscInt                maxChild;
670af33a6ddSJed Brown   CharacteristicPointDA2D temp;
671af33a6ddSJed Brown 
672af33a6ddSJed Brown   PetscFunctionBegin;
673af33a6ddSJed Brown   while ((root * 2 <= bottom) && (!done)) {
674af33a6ddSJed Brown     if (root * 2 == bottom) maxChild = root * 2;
675af33a6ddSJed Brown     else if (queue[root * 2].proc > queue[root * 2 + 1].proc) maxChild = root * 2;
676af33a6ddSJed Brown     else maxChild = root * 2 + 1;
677af33a6ddSJed Brown 
678af33a6ddSJed Brown     if (queue[root].proc < queue[maxChild].proc) {
679af33a6ddSJed Brown       temp            = queue[root];
680af33a6ddSJed Brown       queue[root]     = queue[maxChild];
681af33a6ddSJed Brown       queue[maxChild] = temp;
682af33a6ddSJed Brown       root            = maxChild;
683db4deed7SKarl Rupp     } else done = PETSC_TRUE;
684af33a6ddSJed Brown   }
6853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
686af33a6ddSJed Brown }
687af33a6ddSJed Brown 
688af33a6ddSJed Brown /* [center, left, top-left, top, top-right, right, bottom-right, bottom, bottom-left] */
689*66976f2fSJacob Faibussowitsch static PetscErrorCode DMDAGetNeighborsRank(DM da, PetscMPIInt neighbors[])
690d71ae5a4SJacob Faibussowitsch {
691bff4a2f0SMatthew G. Knepley   DMBoundaryType bx, by;
692af33a6ddSJed Brown   PetscBool      IPeriodic = PETSC_FALSE, JPeriodic = PETSC_FALSE;
693af33a6ddSJed Brown   MPI_Comm       comm;
694af33a6ddSJed Brown   PetscMPIInt    rank;
695af33a6ddSJed Brown   PetscInt     **procs, pi, pj, pim, pip, pjm, pjp, PI, PJ;
696af33a6ddSJed Brown 
697af33a6ddSJed Brown   PetscFunctionBegin;
6989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
6999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7009566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, &PI, &PJ, NULL, NULL, NULL, &bx, &by, NULL, NULL));
701af33a6ddSJed Brown 
702bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC) IPeriodic = PETSC_TRUE;
703bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC) JPeriodic = PETSC_TRUE;
704af33a6ddSJed Brown 
705af33a6ddSJed Brown   neighbors[0] = rank;
706af33a6ddSJed Brown   rank         = 0;
7079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PJ, &procs));
708af33a6ddSJed Brown   for (pj = 0; pj < PJ; pj++) {
7099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PI, &(procs[pj])));
710af33a6ddSJed Brown     for (pi = 0; pi < PI; pi++) {
711af33a6ddSJed Brown       procs[pj][pi] = rank;
712af33a6ddSJed Brown       rank++;
713af33a6ddSJed Brown     }
714af33a6ddSJed Brown   }
715af33a6ddSJed Brown 
716af33a6ddSJed Brown   pi  = neighbors[0] % PI;
717af33a6ddSJed Brown   pj  = neighbors[0] / PI;
7189371c9d4SSatish Balay   pim = pi - 1;
7199371c9d4SSatish Balay   if (pim < 0) pim = PI - 1;
720af33a6ddSJed Brown   pip = (pi + 1) % PI;
7219371c9d4SSatish Balay   pjm = pj - 1;
7229371c9d4SSatish Balay   if (pjm < 0) pjm = PJ - 1;
723af33a6ddSJed Brown   pjp = (pj + 1) % PJ;
724af33a6ddSJed Brown 
725af33a6ddSJed Brown   neighbors[1] = procs[pj][pim];
726af33a6ddSJed Brown   neighbors[2] = procs[pjp][pim];
727af33a6ddSJed Brown   neighbors[3] = procs[pjp][pi];
728af33a6ddSJed Brown   neighbors[4] = procs[pjp][pip];
729af33a6ddSJed Brown   neighbors[5] = procs[pj][pip];
730af33a6ddSJed Brown   neighbors[6] = procs[pjm][pip];
731af33a6ddSJed Brown   neighbors[7] = procs[pjm][pi];
732af33a6ddSJed Brown   neighbors[8] = procs[pjm][pim];
733af33a6ddSJed Brown 
734af33a6ddSJed Brown   if (!IPeriodic) {
735af33a6ddSJed Brown     if (pi == 0) neighbors[1] = neighbors[2] = neighbors[8] = neighbors[0];
736af33a6ddSJed Brown     if (pi == PI - 1) neighbors[4] = neighbors[5] = neighbors[6] = neighbors[0];
737af33a6ddSJed Brown   }
738af33a6ddSJed Brown 
739af33a6ddSJed Brown   if (!JPeriodic) {
740af33a6ddSJed Brown     if (pj == 0) neighbors[6] = neighbors[7] = neighbors[8] = neighbors[0];
741af33a6ddSJed Brown     if (pj == PJ - 1) neighbors[2] = neighbors[3] = neighbors[4] = neighbors[0];
742af33a6ddSJed Brown   }
743af33a6ddSJed Brown 
74448a46eb9SPierre Jolivet   for (pj = 0; pj < PJ; pj++) PetscCall(PetscFree(procs[pj]));
7459566063dSJacob Faibussowitsch   PetscCall(PetscFree(procs));
7463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
747af33a6ddSJed Brown }
748af33a6ddSJed Brown 
749af33a6ddSJed Brown /*
750af33a6ddSJed Brown   SUBDOMAIN NEIGHBORHOOD PROCESS MAP:
751af33a6ddSJed Brown     2 | 3 | 4
752af33a6ddSJed Brown     __|___|__
753af33a6ddSJed Brown     1 | 0 | 5
754af33a6ddSJed Brown     __|___|__
755af33a6ddSJed Brown     8 | 7 | 6
756af33a6ddSJed Brown       |   |
757af33a6ddSJed Brown */
758*66976f2fSJacob Faibussowitsch static PetscInt DMDAGetNeighborRelative(DM da, PetscReal ir, PetscReal jr)
759d71ae5a4SJacob Faibussowitsch {
760af33a6ddSJed Brown   DMDALocalInfo info;
7611cc8b206SBarry Smith   PetscReal     is, ie, js, je;
762af33a6ddSJed Brown 
7639566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
7649371c9d4SSatish Balay   is = (PetscReal)info.xs - 0.5;
7659371c9d4SSatish Balay   ie = (PetscReal)info.xs + info.xm - 0.5;
7669371c9d4SSatish Balay   js = (PetscReal)info.ys - 0.5;
7679371c9d4SSatish Balay   je = (PetscReal)info.ys + info.ym - 0.5;
768af33a6ddSJed Brown 
769af33a6ddSJed Brown   if (ir >= is && ir <= ie) { /* center column */
770bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 0;
771bbd56ea5SKarl Rupp     else if (jr < js) return 7;
772bbd56ea5SKarl Rupp     else return 3;
773af33a6ddSJed Brown   } else if (ir < is) { /* left column */
774bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 1;
775bbd56ea5SKarl Rupp     else if (jr < js) return 8;
776bbd56ea5SKarl Rupp     else return 2;
777af33a6ddSJed Brown   } else { /* right column */
778bbd56ea5SKarl Rupp     if (jr >= js && jr <= je) return 5;
779bbd56ea5SKarl Rupp     else if (jr < js) return 6;
780bbd56ea5SKarl Rupp     else return 4;
781af33a6ddSJed Brown   }
782af33a6ddSJed Brown }
783